xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision ec07b8f88e55a18173750f8ae627228d2c6552c4)
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 
113cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
114cf742fccSHong Zhang {
115cf742fccSHong Zhang   PetscErrorCode    ierr;
116cf742fccSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
117cf742fccSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
118cf742fccSHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
119cf742fccSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
120cf742fccSHong Zhang   Mat               AP_loc,C_loc,C_oth;
1215ca39647SHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz;
122cf742fccSHong Zhang   PetscScalar       *apa;
123cf742fccSHong Zhang   const PetscInt    *cols;
124cf742fccSHong Zhang   const PetscScalar *vals;
125cf742fccSHong Zhang 
126cf742fccSHong Zhang   PetscFunctionBegin;
127cf742fccSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
128cf742fccSHong Zhang 
129cf742fccSHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
130cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
131cf742fccSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
132cf742fccSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
133cf742fccSHong Zhang   }
134cf742fccSHong Zhang 
135cf742fccSHong Zhang   /* 2) get AP_loc */
136cf742fccSHong Zhang   AP_loc = ptap->AP_loc;
137cf742fccSHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
138cf742fccSHong Zhang 
139cf742fccSHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
140cf742fccSHong Zhang   /*-----------------------------------------------------*/
141cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
142cf742fccSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
143cf742fccSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
144cf742fccSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
145cf742fccSHong Zhang   }
146cf742fccSHong Zhang 
147cf742fccSHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
148cf742fccSHong Zhang   /* ---------------------------------------------- */
149cf742fccSHong Zhang   /* get data from symbolic products */
150cf742fccSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
151cf742fccSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
152cf742fccSHong Zhang   api   = ap->i;
153cf742fccSHong Zhang   apj   = ap->j;
154cf742fccSHong Zhang   for (i=0; i<am; i++) {
155cf742fccSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
156cf742fccSHong Zhang     apnz = api[i+1] - api[i];
157b4e8d1b6SHong Zhang     apa = ap->a + api[i];
158b4e8d1b6SHong Zhang     ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
159b4e8d1b6SHong Zhang     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
160cf742fccSHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
161cf742fccSHong Zhang   }
162cf742fccSHong Zhang 
163cf742fccSHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
164cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
165cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
166cf742fccSHong Zhang 
167cf742fccSHong Zhang   C_loc = ptap->C_loc;
168cf742fccSHong Zhang   C_oth = ptap->C_oth;
169cf742fccSHong Zhang 
170cf742fccSHong Zhang   /* add C_loc and Co to to C */
171cf742fccSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
172cf742fccSHong Zhang 
173cf742fccSHong Zhang   /* C_loc -> C */
174cf742fccSHong Zhang   cm    = C_loc->rmap->N;
175cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
176cf742fccSHong Zhang   cols = c_seq->j;
177cf742fccSHong Zhang   vals = c_seq->a;
178cf742fccSHong Zhang   for (i=0; i<cm; i++) {
179cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
180cf742fccSHong Zhang     row = rstart + i;
181cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
182cf742fccSHong Zhang     cols += ncols; vals += ncols;
183cf742fccSHong Zhang   }
184cf742fccSHong Zhang 
185cf742fccSHong Zhang   /* Co -> C, off-processor part */
186cf742fccSHong Zhang   cm = C_oth->rmap->N;
187cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
188cf742fccSHong Zhang   cols = c_seq->j;
189cf742fccSHong Zhang   vals = c_seq->a;
190cf742fccSHong Zhang   for (i=0; i<cm; i++) {
191cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
192cf742fccSHong Zhang     row = p->garray[i];
193cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
194cf742fccSHong Zhang     cols += ncols; vals += ncols;
195cf742fccSHong Zhang   }
196cf742fccSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
197cf742fccSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
198cf742fccSHong Zhang 
199cf742fccSHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
200cf742fccSHong Zhang   PetscFunctionReturn(0);
201cf742fccSHong Zhang }
202cf742fccSHong Zhang 
203e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
2040d3441aeSHong Zhang {
2050d3441aeSHong Zhang   PetscErrorCode      ierr;
2060d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
207681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
2082259aa2eSHong Zhang   MPI_Comm            comm;
2092259aa2eSHong Zhang   PetscMPIInt         size,rank;
21076db11f6SHong Zhang   Mat                 Cmpi,P_loc,P_oth;
21115a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
212d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
213d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
214f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
21515a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
216d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
21715a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
21815a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
21915a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
220445158ffSHong Zhang   PetscLayout         rowmap;
221445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
222445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
223b4e8d1b6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi;
22467a12041SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
22567a12041SHong Zhang   PetscScalar         *apv;
22678d30b94SHong Zhang   PetscTable          ta;
227aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
2288cb82516SHong Zhang   PetscReal           apfill;
229aa690a28SHong Zhang #endif
23067a12041SHong Zhang 
23167a12041SHong Zhang   PetscFunctionBegin;
23267a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
23367a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
23467a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
235ae5f0867Sstefano_zampini 
236ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
237ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
238ae5f0867Sstefano_zampini   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
239ae5f0867Sstefano_zampini 
240e953e47cSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
241e953e47cSHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
242e953e47cSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
243e953e47cSHong Zhang 
244e953e47cSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
24576db11f6SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
246e953e47cSHong Zhang   /* get P_loc by taking all local rows of P */
24776db11f6SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
24876db11f6SHong Zhang 
24976db11f6SHong Zhang   ptap->P_loc = P_loc;
25076db11f6SHong Zhang   ptap->P_oth = P_oth;
251e953e47cSHong Zhang 
252e953e47cSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
253e953e47cSHong Zhang   /* --------------------------------- */
254e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
255e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
256e953e47cSHong Zhang 
257e953e47cSHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
258e953e47cSHong Zhang   /* ----------------------------------------------------------------- */
25976db11f6SHong Zhang   p_loc  = (Mat_SeqAIJ*)P_loc->data;
26076db11f6SHong Zhang   p_oth  = (Mat_SeqAIJ*)P_oth->data;
261e953e47cSHong Zhang 
262e953e47cSHong Zhang   /* create and initialize a linked list */
263e953e47cSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
26476db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
26576db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
266e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
267e953e47cSHong Zhang 
26876db11f6SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
269e953e47cSHong Zhang 
270e953e47cSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
271e953e47cSHong Zhang   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
272e953e47cSHong Zhang   current_space = free_space;
273e953e47cSHong Zhang   nspacedouble  = 0;
274e953e47cSHong Zhang 
275e953e47cSHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
276e953e47cSHong Zhang   api[0] = 0;
277e953e47cSHong Zhang   for (i=0; i<am; i++) {
278e953e47cSHong Zhang     /* diagonal portion: Ad[i,:]*P */
279e953e47cSHong Zhang     ai = ad->i; pi = p_loc->i;
280e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
281e953e47cSHong Zhang     aj  = ad->j + ai[i];
282e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
283e953e47cSHong Zhang       row  = aj[j];
284e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
285e953e47cSHong Zhang       Jptr = p_loc->j + pi[row];
286e953e47cSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
28776db11f6SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
288e953e47cSHong Zhang     }
289e953e47cSHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
290e953e47cSHong Zhang     ai = ao->i; pi = p_oth->i;
291e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
292e953e47cSHong Zhang     aj  = ao->j + ai[i];
293e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
294e953e47cSHong Zhang       row  = aj[j];
295e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
296e953e47cSHong Zhang       Jptr = p_oth->j + pi[row];
29776db11f6SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
298e953e47cSHong Zhang     }
299e953e47cSHong Zhang     apnz     = lnk[0];
300e953e47cSHong Zhang     api[i+1] = api[i] + apnz;
301e953e47cSHong Zhang 
302e953e47cSHong Zhang     /* if free space is not available, double the total space in the list */
303e953e47cSHong Zhang     if (current_space->local_remaining<apnz) {
304e953e47cSHong Zhang       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
305e953e47cSHong Zhang       nspacedouble++;
306e953e47cSHong Zhang     }
307e953e47cSHong Zhang 
308e953e47cSHong Zhang     /* Copy data into free space, then initialize lnk */
30976db11f6SHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
310e953e47cSHong Zhang 
311e953e47cSHong Zhang     current_space->array           += apnz;
312e953e47cSHong Zhang     current_space->local_used      += apnz;
313e953e47cSHong Zhang     current_space->local_remaining -= apnz;
314e953e47cSHong Zhang   }
315e953e47cSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
316e953e47cSHong Zhang   /* destroy list of free space and other temporary array(s) */
317e953e47cSHong Zhang   ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
318e953e47cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
31976db11f6SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
320e953e47cSHong Zhang 
321e953e47cSHong Zhang   /* Create AP_loc for reuse */
322e953e47cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
323e953e47cSHong Zhang 
324e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
325e953e47cSHong Zhang   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
326e953e47cSHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
327e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
328e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
329e953e47cSHong Zhang 
330e953e47cSHong Zhang   if (api[am]) {
331e953e47cSHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
332e953e47cSHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
333e953e47cSHong Zhang   } else {
334e953e47cSHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
335e953e47cSHong Zhang   }
336e953e47cSHong Zhang #endif
337e953e47cSHong Zhang 
338e953e47cSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
339e953e47cSHong Zhang   /* ------------------------------------ */
340e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
341e953e47cSHong Zhang 
342e953e47cSHong Zhang   /* (3) send coj of C_oth to other processors  */
343e953e47cSHong Zhang   /* ------------------------------------------ */
344e953e47cSHong Zhang   /* determine row ownership */
345e953e47cSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
346e953e47cSHong Zhang   rowmap->n  = pn;
347e953e47cSHong Zhang   rowmap->bs = 1;
348e953e47cSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
349e953e47cSHong Zhang   owners = rowmap->range;
350e953e47cSHong Zhang 
351e953e47cSHong Zhang   /* determine the number of messages to send, their lengths */
352e953e47cSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
353e953e47cSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
354e953e47cSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
355e953e47cSHong Zhang 
356e953e47cSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
357e953e47cSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
358e953e47cSHong Zhang   con   = ptap->C_oth->rmap->n;
359e953e47cSHong Zhang   proc  = 0;
360e953e47cSHong Zhang   for (i=0; i<con; i++) {
361e953e47cSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
362e953e47cSHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
363e953e47cSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
364e953e47cSHong Zhang   }
365e953e47cSHong Zhang 
366e953e47cSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
367e953e47cSHong Zhang   owners_co[0] = 0;
368e953e47cSHong Zhang   nsend        = 0;
369e953e47cSHong Zhang   for (proc=0; proc<size; proc++) {
370e953e47cSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
371e953e47cSHong Zhang     if (len_s[proc]) {
372e953e47cSHong Zhang       nsend++;
373e953e47cSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
374e953e47cSHong Zhang       len         += len_si[proc];
375e953e47cSHong Zhang     }
376e953e47cSHong Zhang   }
377e953e47cSHong Zhang 
378e953e47cSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
379e953e47cSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
380e953e47cSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
381e953e47cSHong Zhang 
382e953e47cSHong Zhang   /* post the Irecv and Isend of coj */
383e953e47cSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
384e953e47cSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
385e953e47cSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
386e953e47cSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
387e953e47cSHong Zhang     if (!len_s[proc]) continue;
388e953e47cSHong Zhang     i    = owners_co[proc];
389e953e47cSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
390e953e47cSHong Zhang     k++;
391e953e47cSHong Zhang   }
392e953e47cSHong Zhang 
393e953e47cSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
394e953e47cSHong Zhang   /* ---------------------------------------- */
395e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
396e953e47cSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
397e953e47cSHong Zhang 
398e953e47cSHong Zhang   /* receives coj are complete */
399e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
400e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
401e953e47cSHong Zhang   }
402e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
403e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
404e953e47cSHong Zhang 
405e953e47cSHong Zhang   /* add received column indices into ta to update Crmax */
406e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
407e953e47cSHong Zhang     Jptr = buf_rj[k];
408e953e47cSHong Zhang     for (j=0; j<len_r[k]; j++) {
409e953e47cSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
410e953e47cSHong Zhang     }
411e953e47cSHong Zhang   }
412e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
413e953e47cSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
414e953e47cSHong Zhang 
415e953e47cSHong Zhang   /* (4) send and recv coi */
416e953e47cSHong Zhang   /*-----------------------*/
417e953e47cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
418e953e47cSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
419e953e47cSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
420e953e47cSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
421e953e47cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
422e953e47cSHong Zhang     if (!len_s[proc]) continue;
423e953e47cSHong Zhang     /* form outgoing message for i-structure:
424e953e47cSHong Zhang          buf_si[0]:                 nrows to be sent
425e953e47cSHong Zhang                [1:nrows]:           row index (global)
426e953e47cSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
427e953e47cSHong Zhang     */
428e953e47cSHong Zhang     /*-------------------------------------------*/
429e953e47cSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
430e953e47cSHong Zhang     buf_si_i    = buf_si + nrows+1;
431e953e47cSHong Zhang     buf_si[0]   = nrows;
432e953e47cSHong Zhang     buf_si_i[0] = 0;
433e953e47cSHong Zhang     nrows       = 0;
434e953e47cSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
435e953e47cSHong Zhang       nzi = coi[i+1] - coi[i];
436e953e47cSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
437e953e47cSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
438e953e47cSHong Zhang       nrows++;
439e953e47cSHong Zhang     }
440e953e47cSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
441e953e47cSHong Zhang     k++;
442e953e47cSHong Zhang     buf_si += len_si[proc];
443e953e47cSHong Zhang   }
444e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
445e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
446e953e47cSHong Zhang   }
447e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
448e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
449e953e47cSHong Zhang 
450e953e47cSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
451e953e47cSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
452e953e47cSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
453e953e47cSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
454b4e8d1b6SHong Zhang 
455e953e47cSHong Zhang   /* (5) compute the local portion of Cmpi      */
456e953e47cSHong Zhang   /* ------------------------------------------ */
457e953e47cSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
458e953e47cSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
459e953e47cSHong Zhang   current_space = free_space;
460e953e47cSHong Zhang 
461e953e47cSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
462e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {
463e953e47cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
464e953e47cSHong Zhang     nrows       = *buf_ri_k[k];
465e953e47cSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
466e953e47cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
467e953e47cSHong Zhang   }
468e953e47cSHong Zhang 
469e953e47cSHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
470cf742fccSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
471e953e47cSHong Zhang   for (i=0; i<pn; i++) {
472e953e47cSHong Zhang     /* add C_loc into Cmpi */
473e953e47cSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
474e953e47cSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
475cf742fccSHong Zhang     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
476e953e47cSHong Zhang 
477e953e47cSHong Zhang     /* add received col data into lnk */
478e953e47cSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
479e953e47cSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
480e953e47cSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
481e953e47cSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
482cf742fccSHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
483e953e47cSHong Zhang         nextrow[k]++; nextci[k]++;
484e953e47cSHong Zhang       }
485e953e47cSHong Zhang     }
486e953e47cSHong Zhang     nzi = lnk[0];
487e953e47cSHong Zhang 
488e953e47cSHong Zhang     /* copy data into free space, then initialize lnk */
489cf742fccSHong Zhang     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
490e953e47cSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
491e953e47cSHong Zhang   }
492e953e47cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
493cf742fccSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
494e953e47cSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
495e953e47cSHong Zhang 
496e953e47cSHong Zhang   /* local sizes and preallocation */
497e953e47cSHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
498e953e47cSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
499e953e47cSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
500e953e47cSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
501e953e47cSHong Zhang 
502e953e47cSHong Zhang   /* members in merge */
503e953e47cSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
504e953e47cSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
505e953e47cSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
506e953e47cSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
507e953e47cSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
508e953e47cSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
509e953e47cSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
510e953e47cSHong Zhang 
511e953e47cSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
512e953e47cSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
513e953e47cSHong Zhang   c->ptap         = ptap;
514e953e47cSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
515e953e47cSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
516e953e47cSHong Zhang 
517e953e47cSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
518e953e47cSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
519e953e47cSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
520e953e47cSHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
521e953e47cSHong Zhang   *C                     = Cmpi;
522e953e47cSHong Zhang   PetscFunctionReturn(0);
523e953e47cSHong Zhang }
524e953e47cSHong Zhang 
525e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
526e953e47cSHong Zhang {
527e953e47cSHong Zhang   PetscErrorCode      ierr;
528e953e47cSHong Zhang   Mat_PtAPMPI         *ptap;
529e953e47cSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
530e953e47cSHong Zhang   MPI_Comm            comm;
531e953e47cSHong Zhang   PetscMPIInt         size,rank;
532e953e47cSHong Zhang   Mat                 Cmpi;
533e953e47cSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
534e953e47cSHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
535e953e47cSHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
536e953e47cSHong Zhang   PetscBT             lnkbt;
537e953e47cSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
538e953e47cSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
539e953e47cSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
540e953e47cSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
541e953e47cSHong Zhang   MPI_Request         *swaits,*rwaits;
542e953e47cSHong Zhang   MPI_Status          *sstatus,rstatus;
543e953e47cSHong Zhang   PetscLayout         rowmap;
544e953e47cSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
545e953e47cSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
546e953e47cSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
547*ec07b8f8SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
548e953e47cSHong Zhang   PetscScalar         *apv;
549e953e47cSHong Zhang   PetscTable          ta;
550e953e47cSHong Zhang #if defined(PETSC_HAVE_HYPRE)
551e953e47cSHong Zhang   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
552e953e47cSHong Zhang   PetscInt            nalg = 3;
553e953e47cSHong Zhang #else
554e953e47cSHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
555e953e47cSHong Zhang   PetscInt            nalg = 2;
556e953e47cSHong Zhang #endif
557b4e8d1b6SHong Zhang   PetscInt            alg = 1; /* set default algorithm */
558e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
559e953e47cSHong Zhang   PetscReal           apfill;
560e953e47cSHong Zhang #endif
561e953e47cSHong Zhang #if defined(PTAP_PROFILE)
562e953e47cSHong Zhang   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
563e953e47cSHong Zhang #endif
564b4e8d1b6SHong Zhang   PetscBool           flg;
565e953e47cSHong Zhang 
566e953e47cSHong Zhang   PetscFunctionBegin;
567b4e8d1b6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
568b4e8d1b6SHong Zhang 
569e953e47cSHong Zhang   /* pick an algorithm */
570ae5f0867Sstefano_zampini   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
571b4e8d1b6SHong Zhang   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr);
572ae5f0867Sstefano_zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
573ae5f0867Sstefano_zampini 
574b4e8d1b6SHong Zhang   if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
575b4e8d1b6SHong Zhang     MatInfo     Ainfo,Pinfo;
576b4e8d1b6SHong Zhang     PetscInt    nz_local;
577b4e8d1b6SHong Zhang     PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
578b4e8d1b6SHong Zhang 
579b4e8d1b6SHong Zhang     ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
580b4e8d1b6SHong Zhang     ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
581b4e8d1b6SHong Zhang     nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
582b4e8d1b6SHong Zhang 
583b4e8d1b6SHong Zhang     if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
584b4e8d1b6SHong Zhang     ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
585b4e8d1b6SHong Zhang 
586b4e8d1b6SHong Zhang     if (alg_scalable) {
587b4e8d1b6SHong Zhang       alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
588b4e8d1b6SHong Zhang       ierr = PetscInfo2(P,"Use scalable algorithm, pN %D, fill*nz_allocated %g\n",pN,fill*nz_local);CHKERRQ(ierr);
589b4e8d1b6SHong Zhang     }
590b4e8d1b6SHong Zhang   }
591b4e8d1b6SHong Zhang   /* ierr = PetscInfo2(P,"Use algorithm %D, pN %D\n",alg,pN);CHKERRQ(ierr); */
592b4e8d1b6SHong Zhang 
593e953e47cSHong Zhang   if (alg == 0) {
594e953e47cSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
595cf742fccSHong Zhang     (*C)->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
596e953e47cSHong Zhang     PetscFunctionReturn(0);
597e953e47cSHong Zhang 
598ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
599ae5f0867Sstefano_zampini   } else if (alg == 2) {
600ae5f0867Sstefano_zampini     /* Use boomerAMGBuildCoarseOperator */
601ae5f0867Sstefano_zampini     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
602ae5f0867Sstefano_zampini     PetscFunctionReturn(0);
603ae5f0867Sstefano_zampini #endif
604ae5f0867Sstefano_zampini   }
605ae5f0867Sstefano_zampini 
6068cb82516SHong Zhang #if defined(PTAP_PROFILE)
60780bb4639SHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
6088cb82516SHong Zhang #endif
6098cb82516SHong Zhang 
610e953e47cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
611e953e47cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
612e953e47cSHong Zhang 
613*ec07b8f8SHong Zhang   if (size > 1) {
614*ec07b8f8SHong Zhang     ao = (Mat_SeqAIJ*)(a->B)->data;
615*ec07b8f8SHong Zhang   }
616*ec07b8f8SHong Zhang 
617e953e47cSHong Zhang   /* create symbolic parallel matrix Cmpi */
618e953e47cSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
619e953e47cSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
620e953e47cSHong Zhang 
621e953e47cSHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
622e953e47cSHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
623e953e47cSHong Zhang 
62415a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
62515a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
62615a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
62715a3b8e2SHong Zhang 
62815a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
62915a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
63015a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
63115a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
63215a3b8e2SHong Zhang 
63367a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
63467a12041SHong Zhang   /* --------------------------------- */
635de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
636de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
6378cb82516SHong Zhang #if defined(PTAP_PROFILE)
638445158ffSHong Zhang   ierr = PetscTime(&t11);CHKERRQ(ierr);
6398cb82516SHong Zhang #endif
64015a3b8e2SHong Zhang 
64167a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
64267a12041SHong Zhang   /* ----------------------------------------------------------------- */
64367a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
644*ec07b8f8SHong Zhang   if (ptap->P_oth) {
64567a12041SHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
646*ec07b8f8SHong Zhang   }
647445158ffSHong Zhang 
6489a6dcab7SHong Zhang   /* create and initialize a linked list */
64945d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
6504b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
6514b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
65278d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
653d0e9a2f7SHong 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); */
65478d30b94SHong Zhang 
65578d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
656445158ffSHong Zhang 
6578cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
658*ec07b8f8SHong Zhang   if (ao) {
659f91af8c7SBarry Smith     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
660*ec07b8f8SHong Zhang   } else {
661*ec07b8f8SHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
662*ec07b8f8SHong Zhang   }
663445158ffSHong Zhang   current_space = free_space;
66467a12041SHong Zhang   nspacedouble  = 0;
66567a12041SHong Zhang 
66667a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
66767a12041SHong Zhang   api[0] = 0;
668445158ffSHong Zhang   for (i=0; i<am; i++) {
66967a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
67067a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
67167a12041SHong Zhang     nzi = ai[i+1] - ai[i];
67267a12041SHong Zhang     aj  = ad->j + ai[i];
673445158ffSHong Zhang     for (j=0; j<nzi; j++) {
674445158ffSHong Zhang       row  = aj[j];
67567a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
67667a12041SHong Zhang       Jptr = p_loc->j + pi[row];
677445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
678445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
679445158ffSHong Zhang     }
68067a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
681*ec07b8f8SHong Zhang     if (ao) {
68267a12041SHong Zhang       ai = ao->i; pi = p_oth->i;
68367a12041SHong Zhang       nzi = ai[i+1] - ai[i];
68467a12041SHong Zhang       aj  = ao->j + ai[i];
685445158ffSHong Zhang       for (j=0; j<nzi; j++) {
686445158ffSHong Zhang         row  = aj[j];
68767a12041SHong Zhang         pnz  = pi[row+1] - pi[row];
68867a12041SHong Zhang         Jptr = p_oth->j + pi[row];
689445158ffSHong Zhang         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
690445158ffSHong Zhang       }
691*ec07b8f8SHong Zhang     }
692445158ffSHong Zhang     apnz     = lnk[0];
693445158ffSHong Zhang     api[i+1] = api[i] + apnz;
694445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
695445158ffSHong Zhang 
696445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
697445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
698f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
699445158ffSHong Zhang       nspacedouble++;
700445158ffSHong Zhang     }
701445158ffSHong Zhang 
702445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
703445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
704445158ffSHong Zhang 
705445158ffSHong Zhang     current_space->array           += apnz;
706445158ffSHong Zhang     current_space->local_used      += apnz;
707445158ffSHong Zhang     current_space->local_remaining -= apnz;
708445158ffSHong Zhang   }
709681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
710445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
711445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
712445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
7139a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
7149a6dcab7SHong Zhang 
715aa690a28SHong Zhang   /* Create AP_loc for reuse */
716445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
717aa690a28SHong Zhang 
718aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
719*ec07b8f8SHong Zhang   if (ao) {
720aa690a28SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
721*ec07b8f8SHong Zhang   } else {
722*ec07b8f8SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
723*ec07b8f8SHong Zhang   }
724aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
725aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
726aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
727aa690a28SHong Zhang 
728aa690a28SHong Zhang   if (api[am]) {
729aa690a28SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
730aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
731aa690a28SHong Zhang   } else {
732aa690a28SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
733aa690a28SHong Zhang   }
734aa690a28SHong Zhang #endif
735aa690a28SHong Zhang 
7368cb82516SHong Zhang #if defined(PTAP_PROFILE)
737445158ffSHong Zhang   ierr = PetscTime(&t12);CHKERRQ(ierr);
7388cb82516SHong Zhang #endif
73915a3b8e2SHong Zhang 
740681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
741681d504bSHong Zhang   /* ------------------------------------ */
74267a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
7438cb82516SHong Zhang #if defined(PTAP_PROFILE)
74480bb4639SHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
7458cb82516SHong Zhang #endif
74615a3b8e2SHong Zhang 
74767a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
74867a12041SHong Zhang   /* ------------------------------------------ */
74915a3b8e2SHong Zhang   /* determine row ownership */
750445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
751445158ffSHong Zhang   rowmap->n  = pn;
752445158ffSHong Zhang   rowmap->bs = 1;
753445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
754445158ffSHong Zhang   owners = rowmap->range;
75515a3b8e2SHong Zhang 
75615a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
7578cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
7588cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
75915a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
76015a3b8e2SHong Zhang 
76167a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
76267a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
76367a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
76415a3b8e2SHong Zhang   proc  = 0;
76567a12041SHong Zhang   for (i=0; i<con; i++) {
76615a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
76715a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
76815a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
76915a3b8e2SHong Zhang   }
77015a3b8e2SHong Zhang 
77115a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
77215a3b8e2SHong Zhang   owners_co[0] = 0;
77367a12041SHong Zhang   nsend        = 0;
77415a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
77515a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
77615a3b8e2SHong Zhang     if (len_s[proc]) {
777445158ffSHong Zhang       nsend++;
77815a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
77915a3b8e2SHong Zhang       len         += len_si[proc];
78015a3b8e2SHong Zhang     }
78115a3b8e2SHong Zhang   }
78215a3b8e2SHong Zhang 
78315a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
784445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
785445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
78615a3b8e2SHong Zhang 
78715a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
78815a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
789445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
790445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
79115a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
79215a3b8e2SHong Zhang     if (!len_s[proc]) continue;
79315a3b8e2SHong Zhang     i    = owners_co[proc];
79415a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
79515a3b8e2SHong Zhang     k++;
79615a3b8e2SHong Zhang   }
79715a3b8e2SHong Zhang 
798681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
799681d504bSHong Zhang   /* ---------------------------------------- */
800681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
801681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
802681d504bSHong Zhang 
803681d504bSHong Zhang   /* receives coj are complete */
804445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
805445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
80615a3b8e2SHong Zhang   }
80715a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
808445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
80915a3b8e2SHong Zhang 
81078d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
81178d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
81278d30b94SHong Zhang     Jptr = buf_rj[k];
81378d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
81478d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
81578d30b94SHong Zhang     }
81678d30b94SHong Zhang   }
81778d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
81878d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
8199a6dcab7SHong Zhang 
82015a3b8e2SHong Zhang   /* (4) send and recv coi */
82115a3b8e2SHong Zhang   /*-----------------------*/
82215a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
823445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
82415a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
82515a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
82615a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
82715a3b8e2SHong Zhang     if (!len_s[proc]) continue;
82815a3b8e2SHong Zhang     /* form outgoing message for i-structure:
82915a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
83015a3b8e2SHong Zhang                [1:nrows]:           row index (global)
83115a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
83215a3b8e2SHong Zhang     */
83315a3b8e2SHong Zhang     /*-------------------------------------------*/
83415a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
83515a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
83615a3b8e2SHong Zhang     buf_si[0]   = nrows;
83715a3b8e2SHong Zhang     buf_si_i[0] = 0;
83815a3b8e2SHong Zhang     nrows       = 0;
83915a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
84015a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
84115a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
84215a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
84315a3b8e2SHong Zhang       nrows++;
84415a3b8e2SHong Zhang     }
84515a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
84615a3b8e2SHong Zhang     k++;
84715a3b8e2SHong Zhang     buf_si += len_si[proc];
84815a3b8e2SHong Zhang   }
849681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
850445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
85115a3b8e2SHong Zhang   }
85215a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
853445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
85415a3b8e2SHong Zhang 
8558cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
85615a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
85715a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
85815a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
8598cb82516SHong Zhang #if defined(PTAP_PROFILE)
86080bb4639SHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
8618cb82516SHong Zhang #endif
86267a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
86367a12041SHong Zhang   /* ------------------------------------------ */
86478d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
86578d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
86615a3b8e2SHong Zhang   current_space = free_space;
86715a3b8e2SHong Zhang 
868445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
869445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
87015a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
87115a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
87215a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
87315a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
87415a3b8e2SHong Zhang   }
875445158ffSHong Zhang 
87678d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
87778d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
87815a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
87967a12041SHong Zhang     /* add C_loc into Cmpi */
88067a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
88167a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
88267a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
88315a3b8e2SHong Zhang 
88415a3b8e2SHong Zhang     /* add received col data into lnk */
885445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
88615a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
88715a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
88815a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
88915a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
89015a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
89115a3b8e2SHong Zhang       }
89215a3b8e2SHong Zhang     }
893d0e9a2f7SHong Zhang     nzi = lnk[0];
8948cb82516SHong Zhang 
89515a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
896d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
897d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
89815a3b8e2SHong Zhang   }
89915a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
90015a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
901445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
9028cb82516SHong Zhang #if defined(PTAP_PROFILE)
90380bb4639SHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
9048cb82516SHong Zhang #endif
90580bb4639SHong Zhang 
906ae5f0867Sstefano_zampini   /* local sizes and preallocation */
90715a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
90815a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
90915a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
91015a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
91115a3b8e2SHong Zhang 
912445158ffSHong Zhang   /* members in merge */
913445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
914445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
915445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
916445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
917445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
918445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
919445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
92015a3b8e2SHong Zhang 
92115a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
92215a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
92315a3b8e2SHong Zhang   c->ptap         = ptap;
9241a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
9251a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
92615a3b8e2SHong Zhang 
92778d30b94SHong Zhang   if (alg == 1) {
9288cb82516SHong Zhang     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
92978d30b94SHong Zhang   }
9302259aa2eSHong Zhang 
9311a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9321a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9331a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
934aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
9351a47ec54SHong Zhang   *C                     = Cmpi;
9361a47ec54SHong Zhang 
9378cb82516SHong Zhang #if defined(PTAP_PROFILE)
93880bb4639SHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
939e9c1f85fSHong Zhang   if (rank == 1) {
940445158ffSHong 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);
941e9c1f85fSHong Zhang   }
9428cb82516SHong Zhang #endif
9430d3441aeSHong Zhang   PetscFunctionReturn(0);
9440d3441aeSHong Zhang }
9450d3441aeSHong Zhang 
946aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
9470d3441aeSHong Zhang {
9480d3441aeSHong Zhang   PetscErrorCode    ierr;
9490d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
9500d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
9518cb82516SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
9520d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
9539ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
9540d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
9558cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
9568cb82516SHong Zhang   PetscScalar       *apa;
9570d3441aeSHong Zhang   const PetscInt    *cols;
9580d3441aeSHong Zhang   const PetscScalar *vals;
9598cb82516SHong Zhang #if defined(PTAP_PROFILE)
9608cb82516SHong Zhang   PetscMPIInt       rank;
9618cb82516SHong Zhang   MPI_Comm          comm;
9620d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
9638cb82516SHong Zhang #endif
9640d3441aeSHong Zhang 
9650d3441aeSHong Zhang   PetscFunctionBegin;
9660d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
9670d3441aeSHong Zhang 
968e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
9698cb82516SHong Zhang #if defined(PTAP_PROFILE)
9700d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
9718cb82516SHong Zhang #endif
972748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
9730d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
9740d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
975748c7196SHong Zhang   }
9768cb82516SHong Zhang #if defined(PTAP_PROFILE)
9770d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
9780d3441aeSHong Zhang   eR = t1 - t0;
9798cb82516SHong Zhang #endif
9800d3441aeSHong Zhang 
981e497d3c8SHong Zhang   /* 2) get AP_loc */
9820d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
9838cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
9840d3441aeSHong Zhang 
985e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
9860d3441aeSHong Zhang   /*-----------------------------------------------------*/
987748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
988748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
9890d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
9900d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
991e497d3c8SHong Zhang   }
9920d3441aeSHong Zhang 
9938cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
9948cb82516SHong Zhang   /* ---------------------------------------------- */
9950d3441aeSHong Zhang   /* get data from symbolic products */
9968cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
9978cb82516SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
9988cb82516SHong Zhang   apa   = ptap->apa;
999681d504bSHong Zhang   api   = ap->i;
1000681d504bSHong Zhang   apj   = ap->j;
1001e497d3c8SHong Zhang   for (i=0; i<am; i++) {
1002445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
1003e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
1004e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
1005e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
1006e497d3c8SHong Zhang       col = apj[j+api[i]];
1007e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
1008e497d3c8SHong Zhang       apa[col] = 0.0;
10090d3441aeSHong Zhang     }
1010aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
10110d3441aeSHong Zhang   }
10128cb82516SHong Zhang #if defined(PTAP_PROFILE)
10130d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
10140d3441aeSHong Zhang   eAP = t2 - t1;
10158cb82516SHong Zhang #endif
10160d3441aeSHong Zhang 
10178cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1018445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1019445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
10209ce11a7cSHong Zhang   C_loc = ptap->C_loc;
10219ce11a7cSHong Zhang   C_oth = ptap->C_oth;
10228cb82516SHong Zhang #if defined(PTAP_PROFILE)
10230d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
10240d3441aeSHong Zhang   eCseq = t3 - t2;
10258cb82516SHong Zhang #endif
10260d3441aeSHong Zhang 
10270d3441aeSHong Zhang   /* add C_loc and Co to to C */
10280d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
10290d3441aeSHong Zhang 
10300d3441aeSHong Zhang   /* C_loc -> C */
10310d3441aeSHong Zhang   cm    = C_loc->rmap->N;
10329ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
10338cb82516SHong Zhang   cols = c_seq->j;
10348cb82516SHong Zhang   vals = c_seq->a;
10350d3441aeSHong Zhang   for (i=0; i<cm; i++) {
10369ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
10370d3441aeSHong Zhang     row = rstart + i;
10380d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10398cb82516SHong Zhang     cols += ncols; vals += ncols;
10400d3441aeSHong Zhang   }
10410d3441aeSHong Zhang 
10420d3441aeSHong Zhang   /* Co -> C, off-processor part */
10439ce11a7cSHong Zhang   cm = C_oth->rmap->N;
10449ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
10458cb82516SHong Zhang   cols = c_seq->j;
10468cb82516SHong Zhang   vals = c_seq->a;
10479ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
10489ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
10490d3441aeSHong Zhang     row = p->garray[i];
10500d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10518cb82516SHong Zhang     cols += ncols; vals += ncols;
10520d3441aeSHong Zhang   }
10530d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10540d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10558cb82516SHong Zhang #if defined(PTAP_PROFILE)
10560d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
10570d3441aeSHong Zhang   eCmpi = t4 - t3;
10580d3441aeSHong Zhang 
10598cb82516SHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
10608cb82516SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
10610d3441aeSHong Zhang   if (rank==1) {
10620d3441aeSHong 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);
10630d3441aeSHong Zhang   }
10648cb82516SHong Zhang #endif
1065748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
10660d3441aeSHong Zhang   PetscFunctionReturn(0);
10670d3441aeSHong Zhang }
10680d3441aeSHong Zhang 
10698403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
107042c57489SHong Zhang {
107142c57489SHong Zhang   PetscErrorCode      ierr;
107277584fe6SHong Zhang   Mat                 Cmpi;
1073f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
10740298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1075f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
107642c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
107742c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1078ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
107977584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
1080a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1081d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
108242c57489SHong Zhang   PetscBT             lnkbt;
1083ce94432eSBarry Smith   MPI_Comm            comm;
1084a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
108542c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
108642c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
108724ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
108842c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1089ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
1090ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
109142c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
109277584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
1093d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
109442c57489SHong Zhang 
109542c57489SHong Zhang   PetscFunctionBegin;
1096ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
109783445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
109883445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
109983445d95SHong Zhang 
1100f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1101b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1102b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
11039f4155fbSHong Zhang   ptap->merge = merge;
1104f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
11056c6e00e0SHong Zhang 
11066c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1107b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1108fe615159SHong Zhang 
11096c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
1110a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
11116c6e00e0SHong Zhang 
1112a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1113a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
11146c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
11156c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
111642c57489SHong Zhang 
11172addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
11182addb229SHong Zhang   /*--------------------------------------------------------------------------*/
1119854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
112082412ba7SHong Zhang   api[0] = 0;
112183445d95SHong Zhang 
112283445d95SHong Zhang   /* create and initialize a linked list */
1123a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
112483445d95SHong Zhang 
1125a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
1126f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
1127f4a743e1SHong Zhang   current_space = free_space;
1128f4a743e1SHong Zhang 
1129f4a743e1SHong Zhang   for (i=0; i<am; i++) {
1130f4a743e1SHong Zhang     /* diagonal portion of A */
1131ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
113277584fe6SHong Zhang     aj  = ad->j + adi[i];
1133ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
113477584fe6SHong Zhang       row  = aj[j];
113582412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
1136ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
113783445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
1138a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1139f4a743e1SHong Zhang     }
1140f4a743e1SHong Zhang     /* off-diagonal portion of A */
1141ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
114277584fe6SHong Zhang     aj  = ao->j + aoi[i];
1143ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
114477584fe6SHong Zhang       row  = aj[j];
114582412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
1146ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
1147a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1148f4a743e1SHong Zhang     }
1149a914927fSHong Zhang     apnz     = lnk[0];
115082412ba7SHong Zhang     api[i+1] = api[i] + apnz;
115177584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
1152f4a743e1SHong Zhang 
115383445d95SHong Zhang     /* if free space is not available, double the total space in the list */
115482412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1155f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1156f2b054eeSHong Zhang       nspacedouble++;
1157f4a743e1SHong Zhang     }
1158f4a743e1SHong Zhang 
1159f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
1160a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11612205254eSKarl Rupp 
116282412ba7SHong Zhang     current_space->array           += apnz;
116382412ba7SHong Zhang     current_space->local_used      += apnz;
116482412ba7SHong Zhang     current_space->local_remaining -= apnz;
1165f4a743e1SHong Zhang   }
1166a914927fSHong Zhang 
116782412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
1168f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
1169854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
117077584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
1171118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
1172d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1173f4a743e1SHong Zhang 
11742addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
11752addb229SHong Zhang   /*-----------------------------------------------------------------*/
1176de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
117742c57489SHong Zhang 
1178ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
1179d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
118083445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1181854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1182de0260b3SHong Zhang   coi[0] = 0;
118342c57489SHong Zhang 
1184d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
1185f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am]));
118622d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
118742c57489SHong Zhang   current_space = free_space;
118842c57489SHong Zhang 
1189de0260b3SHong Zhang   for (i=0; i<pon; i++) {
119082412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
119177584fe6SHong Zhang     ptJ = potj + poti[i];
119277584fe6SHong Zhang     for (j=0; j<pnz; j++) {
119377584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
119482412ba7SHong Zhang       apnz = api[row+1] - api[row];
1195ded4f561SHong Zhang       Jptr = apj + api[row];
119683445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1197a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
119842c57489SHong Zhang     }
1199a914927fSHong Zhang     nnz = lnk[0];
120042c57489SHong Zhang 
120183445d95SHong Zhang     /* If free space is not available, double the total space in the list */
1202ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
1203f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1204d16ca5f0SHong Zhang       nspacedouble++;
120542c57489SHong Zhang     }
120642c57489SHong Zhang 
120742c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
1208a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
12092205254eSKarl Rupp 
1210ded4f561SHong Zhang     current_space->array           += nnz;
1211ded4f561SHong Zhang     current_space->local_used      += nnz;
1212ded4f561SHong Zhang     current_space->local_remaining -= nnz;
12132205254eSKarl Rupp 
1214ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
121542c57489SHong Zhang   }
12166b911d16SHong Zhang 
12176b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
1218a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1219118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
1220d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1221de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
122242c57489SHong Zhang 
12236b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
12246b911d16SHong Zhang   /*--------------------------------------------------*/
12256b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
12266b911d16SHong Zhang   len_s        = merge->len_s;
12276b911d16SHong Zhang   merge->nsend = 0;
12286b911d16SHong Zhang 
12296b911d16SHong Zhang 
123042c57489SHong Zhang   /* determine row ownership */
123126283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12327a2fc3feSBarry Smith   merge->rowmap->n  = pn;
12337a2fc3feSBarry Smith   merge->rowmap->bs = 1;
12342205254eSKarl Rupp 
123526283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
12367a2fc3feSBarry Smith   owners = merge->rowmap->range;
123742c57489SHong Zhang 
123842c57489SHong Zhang   /* determine the number of messages to send, their lengths */
1239dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
124083445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1241854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1242de0260b3SHong Zhang 
124383445d95SHong Zhang   proc = 0;
1244de0260b3SHong Zhang   for (i=0; i<pon; i++) {
1245de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
12462addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
1247ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
124883445d95SHong Zhang   }
1249de0260b3SHong Zhang 
12502addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
1251de0260b3SHong Zhang   owners_co[0] = 0;
125242c57489SHong Zhang   for (proc=0; proc<size; proc++) {
1253de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1254ee6e4b5bSHong Zhang     if (len_s[proc]) {
125542c57489SHong Zhang       merge->nsend++;
12562addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
125742c57489SHong Zhang       len         += len_si[proc];
125842c57489SHong Zhang     }
125942c57489SHong Zhang   }
126042c57489SHong Zhang 
1261ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12620298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
126342c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
126442c57489SHong Zhang 
1265ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
1266529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1267529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1268854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
126942c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
127042c57489SHong Zhang     if (!len_s[proc]) continue;
1271de0260b3SHong Zhang     i    = owners_co[proc];
1272529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
127342c57489SHong Zhang     k++;
127442c57489SHong Zhang   }
127542c57489SHong Zhang 
1276ded4f561SHong Zhang   /* receives and sends of coj are complete */
127777584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1278c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1279ded4f561SHong Zhang   }
1280ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
12810c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
128282412ba7SHong Zhang 
12836b911d16SHong Zhang   /* (4) send and recv coi */
12846b911d16SHong Zhang   /*-----------------------*/
1285529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1286529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1287854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
128842c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
128942c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
129042c57489SHong Zhang     if (!len_s[proc]) continue;
129142c57489SHong Zhang     /* form outgoing message for i-structure:
129242c57489SHong Zhang          buf_si[0]:                 nrows to be sent
129342c57489SHong Zhang                [1:nrows]:           row index (global)
129442c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
129542c57489SHong Zhang     */
129642c57489SHong Zhang     /*-------------------------------------------*/
12972addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
129842c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
129942c57489SHong Zhang     buf_si[0]   = nrows;
130042c57489SHong Zhang     buf_si_i[0] = 0;
130142c57489SHong Zhang     nrows       = 0;
1302de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1303ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
1304ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1305de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
130642c57489SHong Zhang       nrows++;
130742c57489SHong Zhang     }
1308529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
130942c57489SHong Zhang     k++;
131042c57489SHong Zhang     buf_si += len_si[proc];
131142c57489SHong Zhang   }
1312ded4f561SHong Zhang   i = merge->nrecv;
1313ded4f561SHong Zhang   while (i--) {
1314c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1315ded4f561SHong Zhang   }
1316ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
13170c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1318a914927fSHong Zhang 
131924ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
132042c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1321ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
132242c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
132342c57489SHong Zhang 
13246b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
13256b911d16SHong Zhang   /*----------------------------------------------*/
1326ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1327ded4f561SHong Zhang 
132824ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
1329854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
133024ecddacSHong Zhang   pti[0] = 0;
133142c57489SHong Zhang 
1332d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1333f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am]));
133422d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
133542c57489SHong Zhang   current_space = free_space;
133642c57489SHong Zhang 
1337dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
133842c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
133942c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
134042c57489SHong Zhang     nrows       = *buf_ri_k[k];
134142c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
134242c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
134342c57489SHong Zhang   }
134442c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
1345936ad1eaSHong Zhang 
134642c57489SHong Zhang   for (i=0; i<pn; i++) {
1347ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
1348ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
134977584fe6SHong Zhang     ptJ = pdtj + pdti[i];
135077584fe6SHong Zhang     for (j=0; j<pnz; j++) {
135177584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1352ded4f561SHong Zhang       apnz = api[row+1] - api[row];
1353ded4f561SHong Zhang       Jptr = apj + api[row];
1354ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1355a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1356ded4f561SHong Zhang     }
1357a914927fSHong Zhang 
135842c57489SHong Zhang     /* add received col data into lnk */
135942c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
136042c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1361ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1362ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1363a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
136442c57489SHong Zhang         nextrow[k]++; nextci[k]++;
136542c57489SHong Zhang       }
136642c57489SHong Zhang     }
1367a914927fSHong Zhang     nnz = lnk[0];
136842c57489SHong Zhang 
136942c57489SHong Zhang     /* if free space is not available, make more free space */
1370ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
1371f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1372d16ca5f0SHong Zhang       nspacedouble++;
137342c57489SHong Zhang     }
137442c57489SHong Zhang     /* copy data into free space, then initialize lnk */
1375a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1376ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13772205254eSKarl Rupp 
1378ded4f561SHong Zhang     current_space->array           += nnz;
1379ded4f561SHong Zhang     current_space->local_used      += nnz;
1380ded4f561SHong Zhang     current_space->local_remaining -= nnz;
13812205254eSKarl Rupp 
138224ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
138342c57489SHong Zhang   }
1384ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
13850572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
138642c57489SHong Zhang 
1387854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
138824ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
138924ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
1390d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
139142c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
139242c57489SHong Zhang 
13936b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
13946b911d16SHong Zhang   /*------------------------------------------*/
139577584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
139677584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
139733d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
139877584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
139977584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
140042c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1401a2f3521dSMark F. Adams 
1402b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
1403b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
1404b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
1405b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
140642c57489SHong Zhang   merge->buf_ri    = buf_ri;
140742c57489SHong Zhang   merge->buf_rj    = buf_rj;
1408de0260b3SHong Zhang   merge->owners_co = owners_co;
140977584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
141042c57489SHong Zhang 
141177584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
141277584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1413f8487c73SHong Zhang   c->ptap     = ptap;
141477584fe6SHong Zhang   ptap->api   = api;
141577584fe6SHong Zhang   ptap->apj   = apj;
14168403a639SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
14178403a639SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
14188403a639SHong Zhang 
14198403a639SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
14208403a639SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
14218403a639SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
14228403a639SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
14238403a639SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
142477584fe6SHong Zhang   *C                     = Cmpi;
1425414894bdSHong Zhang 
1426414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
1427414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
1428414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
1429414894bdSHong Zhang   /* set default scalable */
1430da0a95b2SSatish Balay   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
14312205254eSKarl Rupp 
1432c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
1433414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
14341795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1435414894bdSHong Zhang   } else {
14361795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
1437414894bdSHong Zhang   }
1438414894bdSHong Zhang 
1439f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
144024ecddacSHong Zhang   if (pti[pn] != 0) {
144157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
144257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1443f2b054eeSHong Zhang   } else {
144477584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1445f2b054eeSHong Zhang   }
1446f2b054eeSHong Zhang #endif
144742c57489SHong Zhang   PetscFunctionReturn(0);
144842c57489SHong Zhang }
144942c57489SHong Zhang 
14508403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
145142c57489SHong Zhang {
145242c57489SHong Zhang   PetscErrorCode      ierr;
1453f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
145442c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1455de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
145642c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1457f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
14589f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
14591d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
146082412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
146182412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1462e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1463d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1464ce94432eSBarry Smith   MPI_Comm            comm;
146542c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
146682412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
146742c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
146850f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
146942c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
147042c57489SHong Zhang   MPI_Status          *status;
147182412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
147282412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
1473d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
14746ce94e8fSHong Zhang   PetscBool           scalable;
147538571addSHong Zhang #if defined(PTAP_PROFILE)
1476e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
147738571addSHong Zhang #endif
147842c57489SHong Zhang 
147942c57489SHong Zhang   PetscFunctionBegin;
1480ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
148138571addSHong Zhang #if defined(PTAP_PROFILE)
14828563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
148338571addSHong Zhang #endif
148442c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
148542c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
148642c57489SHong Zhang 
1487f8487c73SHong Zhang   ptap = c->ptap;
1488ce94432eSBarry 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");
1489f8487c73SHong Zhang   merge    = ptap->merge;
1490414894bdSHong Zhang   apa      = ptap->apa;
14916ce94e8fSHong Zhang   scalable = ptap->scalable;
14926c6e00e0SHong Zhang 
1493a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
14946b911d16SHong Zhang   /*-----------------------------------------------------*/
1495f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
14969f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1497f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
14986c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1499b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1500a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
15016c6e00e0SHong Zhang   }
150238571addSHong Zhang #if defined(PTAP_PROFILE)
15038563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
150405d62848SHong Zhang   eP = t1-t0;
150538571addSHong Zhang #endif
150605d62848SHong Zhang   /*
150705d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
150805d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
150905d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
151005d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
151105d62848SHong Zhang    */
1512f8487c73SHong Zhang 
1513f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1514f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
151542c57489SHong Zhang   /* get data from symbolic products */
1516a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1517a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
151889ae1891SBarry Smith   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
151942c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
152042c57489SHong Zhang 
1521de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
15221795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1523de0260b3SHong Zhang 
1524de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
15257a2fc3feSBarry Smith   owners = merge->rowmap->range;
15261795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
152742c57489SHong Zhang 
1528a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1529d9824c15SHong Zhang 
15309516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1531b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
15328cb82516SHong Zhang #if defined(PTAP_PROFILE)
153305d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
15348cb82516SHong Zhang #endif
1535a9d06231SHong Zhang     for (i=0; i<am; i++) {
1536b091e509SHong Zhang #if defined(PTAP_PROFILE)
15378563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1538b091e509SHong Zhang #endif
1539a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1540a9d06231SHong Zhang       /*------------------------------------------------------------*/
1541a9d06231SHong Zhang       apJ = apj + api[i];
1542a9d06231SHong Zhang 
1543a9d06231SHong Zhang       /* diagonal portion of A */
1544a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1545a9d06231SHong Zhang       adj = ad->j + adi[i];
1546a9d06231SHong Zhang       ada = ad->a + adi[i];
1547a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1548a9d06231SHong Zhang         row = adj[j];
1549a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1550a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1551a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1552a9d06231SHong Zhang 
1553a9d06231SHong Zhang         /* perform dense axpy */
1554e900a757SHong Zhang         valtmp = ada[j];
1555a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1556e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1557a9d06231SHong Zhang         }
1558a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1559a9d06231SHong Zhang       }
1560a9d06231SHong Zhang 
1561a9d06231SHong Zhang       /* off-diagonal portion of A */
1562a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1563a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1564a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1565a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1566a9d06231SHong Zhang         row = aoj[j];
1567a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1568a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1569a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1570a9d06231SHong Zhang 
1571a9d06231SHong Zhang         /* perform dense axpy */
1572e900a757SHong Zhang         valtmp = aoa[j];
1573a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1574e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1575a9d06231SHong Zhang         }
1576a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1577a9d06231SHong Zhang       }
1578b091e509SHong Zhang #if defined(PTAP_PROFILE)
15798563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1580b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1581b091e509SHong Zhang #endif
1582a9d06231SHong Zhang 
1583a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1584a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1585a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1586a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1587a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1588e900a757SHong Zhang       poJ = po->j + po->i[i];
1589a9d06231SHong Zhang       pA  = po->a + po->i[i];
1590a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1591e900a757SHong Zhang         row = poJ[j];
1592e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1593e900a757SHong Zhang         cj  = coj + coi[row];
1594e900a757SHong Zhang         ca  = coa + coi[row];
1595a9d06231SHong Zhang         /* perform dense axpy */
1596e900a757SHong Zhang         valtmp = pA[j];
1597a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1598e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1599a9d06231SHong Zhang         }
1600a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1601a9d06231SHong Zhang       }
1602a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1603a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1604e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1605a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1606a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1607e900a757SHong Zhang         row = pdJ[j];
1608e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1609e900a757SHong Zhang         cj  = bj + bi[row];
1610e900a757SHong Zhang         ca  = ba + bi[row];
1611a9d06231SHong Zhang         /* perform dense axpy */
1612e900a757SHong Zhang         valtmp = pA[j];
1613a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1614e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1615a9d06231SHong Zhang         }
1616a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1617a9d06231SHong Zhang       }
16188403a639SHong Zhang 
1619a9d06231SHong Zhang       /* zero the current row of A*P */
1620a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1621b091e509SHong Zhang #if defined(PTAP_PROFILE)
16228563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
162305d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1624b091e509SHong Zhang #endif
1625a9d06231SHong Zhang     }
162638571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1627b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
162838571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1629a9d06231SHong Zhang     pA=pa_loc;
163042c57489SHong Zhang     for (i=0; i<am; i++) {
1631b091e509SHong Zhang #if defined(PTAP_PROFILE)
16328563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1633b091e509SHong Zhang #endif
1634f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
163582412ba7SHong Zhang       apnz = api[i+1] - api[i];
163682412ba7SHong Zhang       apJ  = apj + api[i];
163742c57489SHong Zhang       /* diagonal portion of A */
163882412ba7SHong Zhang       anz = adi[i+1] - adi[i];
16391d633602SHong Zhang       adj = ad->j + adi[i];
16401d633602SHong Zhang       ada = ad->a + adi[i];
164182412ba7SHong Zhang       for (j=0; j<anz; j++) {
16421d633602SHong Zhang         row    = adj[j];
164382412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
164482412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
164582412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1646e900a757SHong Zhang         valtmp = ada[j];
1647f4a743e1SHong Zhang         nextp  = 0;
164882412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
164982412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1650e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
165142c57489SHong Zhang           }
165242c57489SHong Zhang         }
1653dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
165442c57489SHong Zhang       }
165542c57489SHong Zhang       /* off-diagonal portion of A */
165682412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
16571d633602SHong Zhang       aoj = ao->j + aoi[i];
16581d633602SHong Zhang       aoa = ao->a + aoi[i];
165982412ba7SHong Zhang       for (j=0; j<anz; j++) {
16601d633602SHong Zhang         row    = aoj[j];
166182412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
166282412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
166382412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1664e900a757SHong Zhang         valtmp = aoa[j];
1665f4a743e1SHong Zhang         nextp  = 0;
166682412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
166782412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1668e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
166942c57489SHong Zhang           }
167042c57489SHong Zhang         }
1671dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
167242c57489SHong Zhang       }
1673b091e509SHong Zhang #if defined(PTAP_PROFILE)
16748563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1675b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1676b091e509SHong Zhang #endif
167742c57489SHong Zhang 
1678a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1679a9d06231SHong Zhang       /*--------------------------------------------------------------*/
168082412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1681f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
168282412ba7SHong Zhang       for (j=0; j<pnz; j++) {
168342c57489SHong Zhang         nextap = 0;
1684f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
168582412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1686e900a757SHong Zhang           row = *poJ;
1687e900a757SHong Zhang           cj  = coj + coi[row];
1688e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
168982412ba7SHong Zhang         } else {                            /* put the value into Cd */
1690e900a757SHong Zhang           row = *pdJ;
1691e900a757SHong Zhang           cj  = bj + bi[row];
1692e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
169342c57489SHong Zhang         }
1694e900a757SHong Zhang         valtmp = pA[j];
169582412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1696e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
169742c57489SHong Zhang         }
1698dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1699de0260b3SHong Zhang       }
17000496f32cSHong Zhang       pA += pnz;
170142c57489SHong Zhang       /* zero the current row info for A*P */
170282412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1703b091e509SHong Zhang #if defined(PTAP_PROFILE)
17048563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
170505d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1706b091e509SHong Zhang #endif
170742c57489SHong Zhang     }
1708d9824c15SHong Zhang   }
170938571addSHong Zhang #if defined(PTAP_PROFILE)
17108563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
171138571addSHong Zhang #endif
171242c57489SHong Zhang 
1713a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1714a9d06231SHong Zhang   /*------------------------------------*/
171542c57489SHong Zhang   buf_ri = merge->buf_ri;
171642c57489SHong Zhang   buf_rj = merge->buf_rj;
171742c57489SHong Zhang   len_s  = merge->len_s;
1718fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
171942c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
172042c57489SHong Zhang 
1721dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
172242c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
172342c57489SHong Zhang     if (!len_s[proc]) continue;
1724de0260b3SHong Zhang     i    = merge->owners_co[proc];
1725de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
172642c57489SHong Zhang     k++;
172742c57489SHong Zhang   }
17280c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
17290c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
173042c57489SHong Zhang 
1731a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
173242c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
173382412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
173438571addSHong Zhang #if defined(PTAP_PROFILE)
17358563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
173638571addSHong Zhang #endif
173742c57489SHong Zhang 
1738a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1739a9d06231SHong Zhang   /*------------------------------------------------------*/
1740dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
174142c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
174242c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
174342c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
174442c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
174582412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
174642c57489SHong Zhang   }
174742c57489SHong Zhang 
1748de0260b3SHong Zhang   for (i=0; i<cm; i++) {
174982412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
175042c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1751de0260b3SHong Zhang     ba_i = ba + bi[i];
175282412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
175342c57489SHong Zhang     /* add received vals into ba */
175442c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
175542c57489SHong Zhang       /* i-th row */
175642c57489SHong Zhang       if (i == *nextrow[k]) {
175782412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
175882412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
175982412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
176082412ba7SHong Zhang         nextcj = 0;
176182412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
176282412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
176382412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
176442c57489SHong Zhang           }
176542c57489SHong Zhang         }
176682412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1767c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
176842c57489SHong Zhang       }
176942c57489SHong Zhang     }
177082412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
177142c57489SHong Zhang   }
177242c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177342c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177442c57489SHong Zhang 
1775de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1776c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
177742c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
17780572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
177938571addSHong Zhang #if defined(PTAP_PROFILE)
17808563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
17819516a85cSHong Zhang   if (rank==1) {
178205d62848SHong 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);
17839516a85cSHong Zhang   }
178438571addSHong Zhang #endif
178542c57489SHong Zhang   PetscFunctionReturn(0);
178642c57489SHong Zhang }
1787