xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 68eacb73b84ae7f3fd7363217d47f23a8f967155)
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 
15f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
16cc6cb787SHong Zhang {
17cc6cb787SHong Zhang   PetscErrorCode ierr;
18f8487c73SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
19f8487c73SHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
20cc6cb787SHong Zhang 
21cc6cb787SHong Zhang   PetscFunctionBegin;
22f8487c73SHong Zhang   if (ptap) {
23c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI *merge=ptap->merge;
24b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
25f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
26a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
27a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
28c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
2905d62848SHong Zhang     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
3005d62848SHong Zhang     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
318403a639SHong Zhang     if (ptap->AP_loc) { /* used by alg_rap */
32681d504bSHong Zhang       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
33681d504bSHong Zhang       ierr = PetscFree(ap->i);CHKERRQ(ierr);
34681d504bSHong Zhang       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
350d3441aeSHong Zhang       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
368403a639SHong Zhang     } else { /* used by alg_ptap */
378403a639SHong Zhang       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
388403a639SHong Zhang       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
39681d504bSHong Zhang     }
402259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
412259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
42414894bdSHong Zhang     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
43a560ef98SHong Zhang 
448403a639SHong Zhang     if (merge) { /* used by alg_ptap */
45cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
46cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
47cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
48cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
49cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
50c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
51cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
52c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
53cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
54445158ffSHong Zhang       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
55445158ffSHong Zhang       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
5605b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
576bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
58f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
59bf0cc555SLisandro Dalcin     }
60a560ef98SHong Zhang 
61a560ef98SHong Zhang     ierr = ptap->destroy(A);CHKERRQ(ierr);
62f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
63c8b0795cSMark F. Adams   }
64cc6cb787SHong Zhang   PetscFunctionReturn(0);
65cc6cb787SHong Zhang }
66cc6cb787SHong Zhang 
67aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
681a47ec54SHong Zhang {
691a47ec54SHong Zhang   PetscErrorCode ierr;
701a47ec54SHong Zhang   Mat_MPIAIJ     *a     = (Mat_MPIAIJ*)A->data;
711a47ec54SHong Zhang   Mat_PtAPMPI    *ptap  = a->ptap;
721a47ec54SHong Zhang 
731a47ec54SHong Zhang   PetscFunctionBegin;
741a47ec54SHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
751a47ec54SHong Zhang   (*M)->ops->destroy   = ptap->destroy;
761a47ec54SHong Zhang   (*M)->ops->duplicate = ptap->duplicate;
771a47ec54SHong Zhang   PetscFunctionReturn(0);
781a47ec54SHong Zhang }
791a47ec54SHong Zhang 
80150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
8142c57489SHong Zhang {
8242c57489SHong Zhang   PetscErrorCode ierr;
838403a639SHong Zhang   PetscBool      rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */
8467a12041SHong Zhang   MPI_Comm       comm;
8542c57489SHong Zhang 
8642c57489SHong Zhang   PetscFunctionBegin;
8767a12041SHong Zhang   /* check if matrix local sizes are compatible */
8867a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
896c4ed002SBarry 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);
906c4ed002SBarry 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);
9167a12041SHong Zhang 
92c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-matrap",&rap,NULL);CHKERRQ(ierr);
93cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
943ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
958403a639SHong Zhang     if (rap) { /* do R=P^T locally, then C=R*A*P */
96cf3ca8ceSHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
978403a639SHong Zhang     } else {       /* do P^T*A*P */
988403a639SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr);
990d3441aeSHong Zhang     }
1003ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1017a7894deSKris Buschelman   }
1023ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1038403a639SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
1043ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
10542c57489SHong Zhang   PetscFunctionReturn(0);
10642c57489SHong Zhang }
10742c57489SHong Zhang 
108ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
109ae5f0867Sstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
110ae5f0867Sstefano_zampini #endif
111ae5f0867Sstefano_zampini 
112cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
113cf742fccSHong Zhang {
114cf742fccSHong Zhang   PetscErrorCode    ierr;
115cf742fccSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
116cf742fccSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
117cf742fccSHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
118cf742fccSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
119cf742fccSHong Zhang   Mat               AP_loc,C_loc,C_oth;
1205ca39647SHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz;
121cf742fccSHong Zhang   PetscScalar       *apa;
122cf742fccSHong Zhang   const PetscInt    *cols;
123cf742fccSHong Zhang   const PetscScalar *vals;
124cf742fccSHong Zhang 
125cf742fccSHong Zhang   PetscFunctionBegin;
126cf742fccSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
127cf742fccSHong Zhang 
128cf742fccSHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
129cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
130cf742fccSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
131cf742fccSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
132cf742fccSHong Zhang   }
133cf742fccSHong Zhang 
134cf742fccSHong Zhang   /* 2) get AP_loc */
135cf742fccSHong Zhang   AP_loc = ptap->AP_loc;
136cf742fccSHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
137cf742fccSHong Zhang 
138cf742fccSHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
139cf742fccSHong Zhang   /*-----------------------------------------------------*/
140cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
141cf742fccSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
142cf742fccSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
143cf742fccSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
144cf742fccSHong Zhang   }
145cf742fccSHong Zhang 
146cf742fccSHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
147cf742fccSHong Zhang   /* ---------------------------------------------- */
148cf742fccSHong Zhang   /* get data from symbolic products */
149cf742fccSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
15052f7967eSHong Zhang   if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
151c072d3e2SSatish Balay   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!");
15252f7967eSHong Zhang 
153cf742fccSHong Zhang   api   = ap->i;
154cf742fccSHong Zhang   apj   = ap->j;
155cf742fccSHong Zhang   for (i=0; i<am; i++) {
156cf742fccSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
157cf742fccSHong Zhang     apnz = api[i+1] - api[i];
158b4e8d1b6SHong Zhang     apa = ap->a + api[i];
159b4e8d1b6SHong Zhang     ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
160b4e8d1b6SHong Zhang     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
161cf742fccSHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
162cf742fccSHong Zhang   }
163cf742fccSHong Zhang 
164cf742fccSHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
165cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
166cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
167cf742fccSHong Zhang 
168cf742fccSHong Zhang   C_loc = ptap->C_loc;
169cf742fccSHong Zhang   C_oth = ptap->C_oth;
170cf742fccSHong Zhang 
171cf742fccSHong Zhang   /* add C_loc and Co to to C */
172cf742fccSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
173cf742fccSHong Zhang 
174cf742fccSHong Zhang   /* C_loc -> C */
175cf742fccSHong Zhang   cm    = C_loc->rmap->N;
176cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
177cf742fccSHong Zhang   cols = c_seq->j;
178cf742fccSHong Zhang   vals = c_seq->a;
179cf742fccSHong Zhang   for (i=0; i<cm; i++) {
180cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
181cf742fccSHong Zhang     row = rstart + i;
182cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
183cf742fccSHong Zhang     cols += ncols; vals += ncols;
184cf742fccSHong Zhang   }
185cf742fccSHong Zhang 
186cf742fccSHong Zhang   /* Co -> C, off-processor part */
187cf742fccSHong Zhang   cm = C_oth->rmap->N;
188cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
189cf742fccSHong Zhang   cols = c_seq->j;
190cf742fccSHong Zhang   vals = c_seq->a;
191cf742fccSHong Zhang   for (i=0; i<cm; i++) {
192cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
193cf742fccSHong Zhang     row = p->garray[i];
194cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
195cf742fccSHong Zhang     cols += ncols; vals += ncols;
196cf742fccSHong Zhang   }
197cf742fccSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
198cf742fccSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199cf742fccSHong Zhang 
200cf742fccSHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
201cf742fccSHong Zhang   PetscFunctionReturn(0);
202cf742fccSHong Zhang }
203cf742fccSHong Zhang 
204e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
2050d3441aeSHong Zhang {
2060d3441aeSHong Zhang   PetscErrorCode      ierr;
2070d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
208681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
2092259aa2eSHong Zhang   MPI_Comm            comm;
2102259aa2eSHong Zhang   PetscMPIInt         size,rank;
21176db11f6SHong Zhang   Mat                 Cmpi,P_loc,P_oth;
21215a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
213d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
214d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
215f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
21615a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
217d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
21815a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
21915a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
22015a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
221445158ffSHong Zhang   PetscLayout         rowmap;
222445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
223445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
224b4e8d1b6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi;
22552f7967eSHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
22667a12041SHong Zhang   PetscScalar         *apv;
22778d30b94SHong Zhang   PetscTable          ta;
228aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
2298cb82516SHong Zhang   PetscReal           apfill;
230aa690a28SHong Zhang #endif
23167a12041SHong Zhang 
23267a12041SHong Zhang   PetscFunctionBegin;
23367a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
23467a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
23567a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
236ae5f0867Sstefano_zampini 
23752f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
23852f7967eSHong Zhang 
239ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
240ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
241ae5f0867Sstefano_zampini   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
242ae5f0867Sstefano_zampini 
243e953e47cSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
244e953e47cSHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
245e953e47cSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
246e953e47cSHong Zhang 
247e953e47cSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
24876db11f6SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
249e953e47cSHong Zhang   /* get P_loc by taking all local rows of P */
25076db11f6SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
25176db11f6SHong Zhang 
25276db11f6SHong Zhang   ptap->P_loc = P_loc;
25376db11f6SHong Zhang   ptap->P_oth = P_oth;
254e953e47cSHong Zhang 
255e953e47cSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
256e953e47cSHong Zhang   /* --------------------------------- */
257e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
258e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
259e953e47cSHong Zhang 
260e953e47cSHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
261e953e47cSHong Zhang   /* ----------------------------------------------------------------- */
26276db11f6SHong Zhang   p_loc  = (Mat_SeqAIJ*)P_loc->data;
26352f7967eSHong Zhang   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
264e953e47cSHong Zhang 
265e953e47cSHong Zhang   /* create and initialize a linked list */
266e953e47cSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
26776db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
26876db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
269e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
270e953e47cSHong Zhang 
27176db11f6SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
272e953e47cSHong Zhang 
273e953e47cSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
27452f7967eSHong Zhang   if (ao) {
275e953e47cSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
27652f7967eSHong Zhang   } else {
27752f7967eSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
27852f7967eSHong Zhang   }
279e953e47cSHong Zhang   current_space = free_space;
280e953e47cSHong Zhang   nspacedouble  = 0;
281e953e47cSHong Zhang 
282e953e47cSHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
283e953e47cSHong Zhang   api[0] = 0;
284e953e47cSHong Zhang   for (i=0; i<am; i++) {
285e953e47cSHong Zhang     /* diagonal portion: Ad[i,:]*P */
286e953e47cSHong Zhang     ai = ad->i; pi = p_loc->i;
287e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
288e953e47cSHong Zhang     aj  = ad->j + ai[i];
289e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
290e953e47cSHong Zhang       row  = aj[j];
291e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
292e953e47cSHong Zhang       Jptr = p_loc->j + pi[row];
293e953e47cSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
29476db11f6SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
295e953e47cSHong Zhang     }
296e953e47cSHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
29752f7967eSHong Zhang     if (ao) {
298e953e47cSHong Zhang       ai = ao->i; pi = p_oth->i;
299e953e47cSHong Zhang       nzi = ai[i+1] - ai[i];
300e953e47cSHong Zhang       aj  = ao->j + ai[i];
301e953e47cSHong Zhang       for (j=0; j<nzi; j++) {
302e953e47cSHong Zhang         row  = aj[j];
303e953e47cSHong Zhang         pnz  = pi[row+1] - pi[row];
304e953e47cSHong Zhang         Jptr = p_oth->j + pi[row];
30576db11f6SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
306e953e47cSHong Zhang       }
30752f7967eSHong Zhang     }
308e953e47cSHong Zhang     apnz     = lnk[0];
309e953e47cSHong Zhang     api[i+1] = api[i] + apnz;
310e953e47cSHong Zhang 
311e953e47cSHong Zhang     /* if free space is not available, double the total space in the list */
312e953e47cSHong Zhang     if (current_space->local_remaining<apnz) {
313e953e47cSHong Zhang       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
314e953e47cSHong Zhang       nspacedouble++;
315e953e47cSHong Zhang     }
316e953e47cSHong Zhang 
317e953e47cSHong Zhang     /* Copy data into free space, then initialize lnk */
31876db11f6SHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
319e953e47cSHong Zhang 
320e953e47cSHong Zhang     current_space->array           += apnz;
321e953e47cSHong Zhang     current_space->local_used      += apnz;
322e953e47cSHong Zhang     current_space->local_remaining -= apnz;
323e953e47cSHong Zhang   }
324e953e47cSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
325e953e47cSHong Zhang   /* destroy list of free space and other temporary array(s) */
326e953e47cSHong Zhang   ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
327e953e47cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
32876db11f6SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
329e953e47cSHong Zhang 
330e953e47cSHong Zhang   /* Create AP_loc for reuse */
331e953e47cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
332e953e47cSHong Zhang 
333e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
33452f7967eSHong Zhang   if (ao) {
335e953e47cSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
33652f7967eSHong Zhang   } else {
33752f7967eSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
33852f7967eSHong Zhang   }
339e953e47cSHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
340e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
341e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
342e953e47cSHong Zhang 
343e953e47cSHong Zhang   if (api[am]) {
344e953e47cSHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
345e953e47cSHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
346e953e47cSHong Zhang   } else {
347e953e47cSHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
348e953e47cSHong Zhang   }
349e953e47cSHong Zhang #endif
350e953e47cSHong Zhang 
351e953e47cSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
352e953e47cSHong Zhang   /* ------------------------------------ */
353e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
354e953e47cSHong Zhang 
355e953e47cSHong Zhang   /* (3) send coj of C_oth to other processors  */
356e953e47cSHong Zhang   /* ------------------------------------------ */
357e953e47cSHong Zhang   /* determine row ownership */
358e953e47cSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
359e953e47cSHong Zhang   rowmap->n  = pn;
360e953e47cSHong Zhang   rowmap->bs = 1;
361e953e47cSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
362e953e47cSHong Zhang   owners = rowmap->range;
363e953e47cSHong Zhang 
364e953e47cSHong Zhang   /* determine the number of messages to send, their lengths */
365e953e47cSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
366e953e47cSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
367e953e47cSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
368e953e47cSHong Zhang 
369e953e47cSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
370e953e47cSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
371e953e47cSHong Zhang   con   = ptap->C_oth->rmap->n;
372e953e47cSHong Zhang   proc  = 0;
373e953e47cSHong Zhang   for (i=0; i<con; i++) {
374e953e47cSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
375e953e47cSHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
376e953e47cSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
377e953e47cSHong Zhang   }
378e953e47cSHong Zhang 
379e953e47cSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
380e953e47cSHong Zhang   owners_co[0] = 0;
381e953e47cSHong Zhang   nsend        = 0;
382e953e47cSHong Zhang   for (proc=0; proc<size; proc++) {
383e953e47cSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
384e953e47cSHong Zhang     if (len_s[proc]) {
385e953e47cSHong Zhang       nsend++;
386e953e47cSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
387e953e47cSHong Zhang       len         += len_si[proc];
388e953e47cSHong Zhang     }
389e953e47cSHong Zhang   }
390e953e47cSHong Zhang 
391e953e47cSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
392e953e47cSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
393e953e47cSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
394e953e47cSHong Zhang 
395e953e47cSHong Zhang   /* post the Irecv and Isend of coj */
396e953e47cSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
397e953e47cSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
398e953e47cSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
399e953e47cSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
400e953e47cSHong Zhang     if (!len_s[proc]) continue;
401e953e47cSHong Zhang     i    = owners_co[proc];
402e953e47cSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
403e953e47cSHong Zhang     k++;
404e953e47cSHong Zhang   }
405e953e47cSHong Zhang 
406e953e47cSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
407e953e47cSHong Zhang   /* ---------------------------------------- */
408e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
409e953e47cSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
410e953e47cSHong Zhang 
411e953e47cSHong Zhang   /* receives coj are complete */
412e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
413e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
414e953e47cSHong Zhang   }
415e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
416e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
417e953e47cSHong Zhang 
418e953e47cSHong Zhang   /* add received column indices into ta to update Crmax */
419e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
420e953e47cSHong Zhang     Jptr = buf_rj[k];
421e953e47cSHong Zhang     for (j=0; j<len_r[k]; j++) {
422e953e47cSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
423e953e47cSHong Zhang     }
424e953e47cSHong Zhang   }
425e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
426e953e47cSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
427e953e47cSHong Zhang 
428e953e47cSHong Zhang   /* (4) send and recv coi */
429e953e47cSHong Zhang   /*-----------------------*/
430e953e47cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
431e953e47cSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
432e953e47cSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
433e953e47cSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
434e953e47cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
435e953e47cSHong Zhang     if (!len_s[proc]) continue;
436e953e47cSHong Zhang     /* form outgoing message for i-structure:
437e953e47cSHong Zhang          buf_si[0]:                 nrows to be sent
438e953e47cSHong Zhang                [1:nrows]:           row index (global)
439e953e47cSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
440e953e47cSHong Zhang     */
441e953e47cSHong Zhang     /*-------------------------------------------*/
442e953e47cSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
443e953e47cSHong Zhang     buf_si_i    = buf_si + nrows+1;
444e953e47cSHong Zhang     buf_si[0]   = nrows;
445e953e47cSHong Zhang     buf_si_i[0] = 0;
446e953e47cSHong Zhang     nrows       = 0;
447e953e47cSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
448e953e47cSHong Zhang       nzi = coi[i+1] - coi[i];
449e953e47cSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
450e953e47cSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
451e953e47cSHong Zhang       nrows++;
452e953e47cSHong Zhang     }
453e953e47cSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
454e953e47cSHong Zhang     k++;
455e953e47cSHong Zhang     buf_si += len_si[proc];
456e953e47cSHong Zhang   }
457e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
458e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
459e953e47cSHong Zhang   }
460e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
461e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
462e953e47cSHong Zhang 
463e953e47cSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
464e953e47cSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
465e953e47cSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
466e953e47cSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
467b4e8d1b6SHong Zhang 
468e953e47cSHong Zhang   /* (5) compute the local portion of Cmpi      */
469e953e47cSHong Zhang   /* ------------------------------------------ */
470e953e47cSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
471e953e47cSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
472e953e47cSHong Zhang   current_space = free_space;
473e953e47cSHong Zhang 
474e953e47cSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
475e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {
476e953e47cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
477e953e47cSHong Zhang     nrows       = *buf_ri_k[k];
478e953e47cSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
479e953e47cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
480e953e47cSHong Zhang   }
481e953e47cSHong Zhang 
482e953e47cSHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
483cf742fccSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
484e953e47cSHong Zhang   for (i=0; i<pn; i++) {
485e953e47cSHong Zhang     /* add C_loc into Cmpi */
486e953e47cSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
487e953e47cSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
488cf742fccSHong Zhang     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
489e953e47cSHong Zhang 
490e953e47cSHong Zhang     /* add received col data into lnk */
491e953e47cSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
492e953e47cSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
493e953e47cSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
494e953e47cSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
495cf742fccSHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
496e953e47cSHong Zhang         nextrow[k]++; nextci[k]++;
497e953e47cSHong Zhang       }
498e953e47cSHong Zhang     }
499e953e47cSHong Zhang     nzi = lnk[0];
500e953e47cSHong Zhang 
501e953e47cSHong Zhang     /* copy data into free space, then initialize lnk */
502cf742fccSHong Zhang     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
503e953e47cSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
504e953e47cSHong Zhang   }
505e953e47cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
506cf742fccSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
507e953e47cSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
508e953e47cSHong Zhang 
509e953e47cSHong Zhang   /* local sizes and preallocation */
510e953e47cSHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
511e953e47cSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
512e953e47cSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
513e953e47cSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
514e953e47cSHong Zhang 
515e953e47cSHong Zhang   /* members in merge */
516e953e47cSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
517e953e47cSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
518e953e47cSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
519e953e47cSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
520e953e47cSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
521e953e47cSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
522e953e47cSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
523e953e47cSHong Zhang 
524e953e47cSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
525e953e47cSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
526e953e47cSHong Zhang   c->ptap         = ptap;
527e953e47cSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
528e953e47cSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
529e953e47cSHong Zhang 
530e953e47cSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
531e953e47cSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
532e953e47cSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
533e953e47cSHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
534e953e47cSHong Zhang   *C                     = Cmpi;
535e953e47cSHong Zhang   PetscFunctionReturn(0);
536e953e47cSHong Zhang }
537e953e47cSHong Zhang 
538e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
539e953e47cSHong Zhang {
540e953e47cSHong Zhang   PetscErrorCode      ierr;
541e953e47cSHong Zhang   Mat_PtAPMPI         *ptap;
542e953e47cSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
543e953e47cSHong Zhang   MPI_Comm            comm;
544e953e47cSHong Zhang   PetscMPIInt         size,rank;
545e953e47cSHong Zhang   Mat                 Cmpi;
546e953e47cSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
547e953e47cSHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
548e953e47cSHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
549e953e47cSHong Zhang   PetscBT             lnkbt;
550e953e47cSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
551e953e47cSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
552e953e47cSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
553e953e47cSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
554e953e47cSHong Zhang   MPI_Request         *swaits,*rwaits;
555e953e47cSHong Zhang   MPI_Status          *sstatus,rstatus;
556e953e47cSHong Zhang   PetscLayout         rowmap;
557e953e47cSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
558e953e47cSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
559e953e47cSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
560ec07b8f8SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
561e953e47cSHong Zhang   PetscScalar         *apv;
562e953e47cSHong Zhang   PetscTable          ta;
563e953e47cSHong Zhang #if defined(PETSC_HAVE_HYPRE)
564e953e47cSHong Zhang   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
565e953e47cSHong Zhang   PetscInt            nalg = 3;
566e953e47cSHong Zhang #else
567e953e47cSHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
568e953e47cSHong Zhang   PetscInt            nalg = 2;
569e953e47cSHong Zhang #endif
570b4e8d1b6SHong Zhang   PetscInt            alg = 1; /* set default algorithm */
571e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
572e953e47cSHong Zhang   PetscReal           apfill;
573e953e47cSHong Zhang #endif
574e953e47cSHong Zhang #if defined(PTAP_PROFILE)
575e953e47cSHong Zhang   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
576e953e47cSHong Zhang #endif
577b4e8d1b6SHong Zhang   PetscBool           flg;
578e953e47cSHong Zhang 
579e953e47cSHong Zhang   PetscFunctionBegin;
580b4e8d1b6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
581b4e8d1b6SHong Zhang 
582e953e47cSHong Zhang   /* pick an algorithm */
583ae5f0867Sstefano_zampini   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
584*68eacb73SHong Zhang   PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
585b4e8d1b6SHong Zhang   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr);
586ae5f0867Sstefano_zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
587ae5f0867Sstefano_zampini 
588b4e8d1b6SHong Zhang   if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
589b4e8d1b6SHong Zhang     MatInfo     Ainfo,Pinfo;
590b4e8d1b6SHong Zhang     PetscInt    nz_local;
591b4e8d1b6SHong Zhang     PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
592b4e8d1b6SHong Zhang 
593b4e8d1b6SHong Zhang     ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
594b4e8d1b6SHong Zhang     ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
595b4e8d1b6SHong Zhang     nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
596b4e8d1b6SHong Zhang 
597b4e8d1b6SHong Zhang     if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
598b4e8d1b6SHong Zhang     ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
599b4e8d1b6SHong Zhang 
600b4e8d1b6SHong Zhang     if (alg_scalable) {
601b4e8d1b6SHong Zhang       alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
602b4e8d1b6SHong Zhang       ierr = PetscInfo2(P,"Use scalable algorithm, pN %D, fill*nz_allocated %g\n",pN,fill*nz_local);CHKERRQ(ierr);
603b4e8d1b6SHong Zhang     }
604b4e8d1b6SHong Zhang   }
605b4e8d1b6SHong Zhang   /* ierr = PetscInfo2(P,"Use algorithm %D, pN %D\n",alg,pN);CHKERRQ(ierr); */
606b4e8d1b6SHong Zhang 
607e953e47cSHong Zhang   if (alg == 0) {
608e953e47cSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
609cf742fccSHong Zhang     (*C)->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
610e953e47cSHong Zhang     PetscFunctionReturn(0);
611e953e47cSHong Zhang 
612ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
613ae5f0867Sstefano_zampini   } else if (alg == 2) {
614ae5f0867Sstefano_zampini     /* Use boomerAMGBuildCoarseOperator */
615ae5f0867Sstefano_zampini     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
616ae5f0867Sstefano_zampini     PetscFunctionReturn(0);
617ae5f0867Sstefano_zampini #endif
618ae5f0867Sstefano_zampini   }
619ae5f0867Sstefano_zampini 
6208cb82516SHong Zhang #if defined(PTAP_PROFILE)
62180bb4639SHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
6228cb82516SHong Zhang #endif
6238cb82516SHong Zhang 
624e953e47cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
625e953e47cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
626e953e47cSHong Zhang 
62752f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
628ec07b8f8SHong Zhang 
629e953e47cSHong Zhang   /* create symbolic parallel matrix Cmpi */
630e953e47cSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
631e953e47cSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
632e953e47cSHong Zhang 
633e953e47cSHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
634e953e47cSHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
635e953e47cSHong Zhang 
63615a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
63715a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
63815a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
63915a3b8e2SHong Zhang 
64015a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
64115a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
64215a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
64315a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
64415a3b8e2SHong Zhang 
64567a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
64667a12041SHong Zhang   /* --------------------------------- */
647de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
648de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
6498cb82516SHong Zhang #if defined(PTAP_PROFILE)
650445158ffSHong Zhang   ierr = PetscTime(&t11);CHKERRQ(ierr);
6518cb82516SHong Zhang #endif
65215a3b8e2SHong Zhang 
65367a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
65467a12041SHong Zhang   /* ----------------------------------------------------------------- */
65567a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
65652f7967eSHong Zhang   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
657445158ffSHong Zhang 
6589a6dcab7SHong Zhang   /* create and initialize a linked list */
65945d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
6604b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
6614b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
66278d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
663d0e9a2f7SHong 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); */
66478d30b94SHong Zhang 
66578d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
666445158ffSHong Zhang 
6678cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
668ec07b8f8SHong Zhang   if (ao) {
669f91af8c7SBarry Smith     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
670ec07b8f8SHong Zhang   } else {
671ec07b8f8SHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
672ec07b8f8SHong Zhang   }
673445158ffSHong Zhang   current_space = free_space;
67467a12041SHong Zhang   nspacedouble  = 0;
67567a12041SHong Zhang 
67667a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
67767a12041SHong Zhang   api[0] = 0;
678445158ffSHong Zhang   for (i=0; i<am; i++) {
67967a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
68067a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
68167a12041SHong Zhang     nzi = ai[i+1] - ai[i];
68267a12041SHong Zhang     aj  = ad->j + ai[i];
683445158ffSHong Zhang     for (j=0; j<nzi; j++) {
684445158ffSHong Zhang       row  = aj[j];
68567a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
68667a12041SHong Zhang       Jptr = p_loc->j + pi[row];
687445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
688445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
689445158ffSHong Zhang     }
69067a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
691ec07b8f8SHong Zhang     if (ao) {
69267a12041SHong Zhang       ai = ao->i; pi = p_oth->i;
69367a12041SHong Zhang       nzi = ai[i+1] - ai[i];
69467a12041SHong Zhang       aj  = ao->j + ai[i];
695445158ffSHong Zhang       for (j=0; j<nzi; j++) {
696445158ffSHong Zhang         row  = aj[j];
69767a12041SHong Zhang         pnz  = pi[row+1] - pi[row];
69867a12041SHong Zhang         Jptr = p_oth->j + pi[row];
699445158ffSHong Zhang         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
700445158ffSHong Zhang       }
701ec07b8f8SHong Zhang     }
702445158ffSHong Zhang     apnz     = lnk[0];
703445158ffSHong Zhang     api[i+1] = api[i] + apnz;
704445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
705445158ffSHong Zhang 
706445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
707445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
708f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
709445158ffSHong Zhang       nspacedouble++;
710445158ffSHong Zhang     }
711445158ffSHong Zhang 
712445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
713445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
714445158ffSHong Zhang 
715445158ffSHong Zhang     current_space->array           += apnz;
716445158ffSHong Zhang     current_space->local_used      += apnz;
717445158ffSHong Zhang     current_space->local_remaining -= apnz;
718445158ffSHong Zhang   }
719681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
720445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
721445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
722445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
7239a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
7249a6dcab7SHong Zhang 
725aa690a28SHong Zhang   /* Create AP_loc for reuse */
726445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
727aa690a28SHong Zhang 
728aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
729ec07b8f8SHong Zhang   if (ao) {
730aa690a28SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
731ec07b8f8SHong Zhang   } else {
732ec07b8f8SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
733ec07b8f8SHong Zhang   }
734aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
735aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
736aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
737aa690a28SHong Zhang 
738aa690a28SHong Zhang   if (api[am]) {
739aa690a28SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
740aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
741aa690a28SHong Zhang   } else {
742aa690a28SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
743aa690a28SHong Zhang   }
744aa690a28SHong Zhang #endif
745aa690a28SHong Zhang 
7468cb82516SHong Zhang #if defined(PTAP_PROFILE)
747445158ffSHong Zhang   ierr = PetscTime(&t12);CHKERRQ(ierr);
7488cb82516SHong Zhang #endif
74915a3b8e2SHong Zhang 
750681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
751681d504bSHong Zhang   /* ------------------------------------ */
75267a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
7538cb82516SHong Zhang #if defined(PTAP_PROFILE)
75480bb4639SHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
7558cb82516SHong Zhang #endif
75615a3b8e2SHong Zhang 
75767a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
75867a12041SHong Zhang   /* ------------------------------------------ */
75915a3b8e2SHong Zhang   /* determine row ownership */
760445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
761445158ffSHong Zhang   rowmap->n  = pn;
762445158ffSHong Zhang   rowmap->bs = 1;
763445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
764445158ffSHong Zhang   owners = rowmap->range;
76515a3b8e2SHong Zhang 
76615a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
7678cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
7688cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
76915a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
77015a3b8e2SHong Zhang 
77167a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
77267a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
77367a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
77415a3b8e2SHong Zhang   proc  = 0;
77567a12041SHong Zhang   for (i=0; i<con; i++) {
77615a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
77715a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
77815a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
77915a3b8e2SHong Zhang   }
78015a3b8e2SHong Zhang 
78115a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
78215a3b8e2SHong Zhang   owners_co[0] = 0;
78367a12041SHong Zhang   nsend        = 0;
78415a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
78515a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
78615a3b8e2SHong Zhang     if (len_s[proc]) {
787445158ffSHong Zhang       nsend++;
78815a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
78915a3b8e2SHong Zhang       len         += len_si[proc];
79015a3b8e2SHong Zhang     }
79115a3b8e2SHong Zhang   }
79215a3b8e2SHong Zhang 
79315a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
794445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
795445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
79615a3b8e2SHong Zhang 
79715a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
79815a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
799445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
800445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
80115a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
80215a3b8e2SHong Zhang     if (!len_s[proc]) continue;
80315a3b8e2SHong Zhang     i    = owners_co[proc];
80415a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
80515a3b8e2SHong Zhang     k++;
80615a3b8e2SHong Zhang   }
80715a3b8e2SHong Zhang 
808681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
809681d504bSHong Zhang   /* ---------------------------------------- */
810681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
811681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
812681d504bSHong Zhang 
813681d504bSHong Zhang   /* receives coj are complete */
814445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
815445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
81615a3b8e2SHong Zhang   }
81715a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
818445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
81915a3b8e2SHong Zhang 
82078d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
82178d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
82278d30b94SHong Zhang     Jptr = buf_rj[k];
82378d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
82478d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
82578d30b94SHong Zhang     }
82678d30b94SHong Zhang   }
82778d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
82878d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
8299a6dcab7SHong Zhang 
83015a3b8e2SHong Zhang   /* (4) send and recv coi */
83115a3b8e2SHong Zhang   /*-----------------------*/
83215a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
833445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
83415a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
83515a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
83615a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
83715a3b8e2SHong Zhang     if (!len_s[proc]) continue;
83815a3b8e2SHong Zhang     /* form outgoing message for i-structure:
83915a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
84015a3b8e2SHong Zhang                [1:nrows]:           row index (global)
84115a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
84215a3b8e2SHong Zhang     */
84315a3b8e2SHong Zhang     /*-------------------------------------------*/
84415a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
84515a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
84615a3b8e2SHong Zhang     buf_si[0]   = nrows;
84715a3b8e2SHong Zhang     buf_si_i[0] = 0;
84815a3b8e2SHong Zhang     nrows       = 0;
84915a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
85015a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
85115a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
85215a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
85315a3b8e2SHong Zhang       nrows++;
85415a3b8e2SHong Zhang     }
85515a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
85615a3b8e2SHong Zhang     k++;
85715a3b8e2SHong Zhang     buf_si += len_si[proc];
85815a3b8e2SHong Zhang   }
859681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
860445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
86115a3b8e2SHong Zhang   }
86215a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
863445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
86415a3b8e2SHong Zhang 
8658cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
86615a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
86715a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
86815a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
8698cb82516SHong Zhang #if defined(PTAP_PROFILE)
87080bb4639SHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
8718cb82516SHong Zhang #endif
87267a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
87367a12041SHong Zhang   /* ------------------------------------------ */
87478d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
87578d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
87615a3b8e2SHong Zhang   current_space = free_space;
87715a3b8e2SHong Zhang 
878445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
879445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
88015a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
88115a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
88215a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
88315a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
88415a3b8e2SHong Zhang   }
885445158ffSHong Zhang 
88678d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
88778d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
88815a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
88967a12041SHong Zhang     /* add C_loc into Cmpi */
89067a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
89167a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
89267a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
89315a3b8e2SHong Zhang 
89415a3b8e2SHong Zhang     /* add received col data into lnk */
895445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
89615a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
89715a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
89815a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
89915a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
90015a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
90115a3b8e2SHong Zhang       }
90215a3b8e2SHong Zhang     }
903d0e9a2f7SHong Zhang     nzi = lnk[0];
9048cb82516SHong Zhang 
90515a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
906d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
907d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
90815a3b8e2SHong Zhang   }
90915a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
91015a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
911445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
9128cb82516SHong Zhang #if defined(PTAP_PROFILE)
91380bb4639SHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
9148cb82516SHong Zhang #endif
91580bb4639SHong Zhang 
916ae5f0867Sstefano_zampini   /* local sizes and preallocation */
91715a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
91815a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
91915a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
92015a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
92115a3b8e2SHong Zhang 
922445158ffSHong Zhang   /* members in merge */
923445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
924445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
925445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
926445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
927445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
928445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
929445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
93015a3b8e2SHong Zhang 
93115a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
93215a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
93315a3b8e2SHong Zhang   c->ptap         = ptap;
9341a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
9351a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
93615a3b8e2SHong Zhang 
93778d30b94SHong Zhang   if (alg == 1) {
9388cb82516SHong Zhang     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
93978d30b94SHong Zhang   }
9402259aa2eSHong Zhang 
9411a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9421a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9431a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
944aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
9451a47ec54SHong Zhang   *C                     = Cmpi;
9461a47ec54SHong Zhang 
9478cb82516SHong Zhang #if defined(PTAP_PROFILE)
94880bb4639SHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
949e9c1f85fSHong Zhang   if (rank == 1) {
950445158ffSHong 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);
951e9c1f85fSHong Zhang   }
9528cb82516SHong Zhang #endif
9530d3441aeSHong Zhang   PetscFunctionReturn(0);
9540d3441aeSHong Zhang }
9550d3441aeSHong Zhang 
956aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
9570d3441aeSHong Zhang {
9580d3441aeSHong Zhang   PetscErrorCode    ierr;
9590d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
9600d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
9612dd9e643SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
9620d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
9639ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
9640d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
9658cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
9668cb82516SHong Zhang   PetscScalar       *apa;
9670d3441aeSHong Zhang   const PetscInt    *cols;
9680d3441aeSHong Zhang   const PetscScalar *vals;
9698cb82516SHong Zhang #if defined(PTAP_PROFILE)
9708cb82516SHong Zhang   PetscMPIInt       rank;
9718cb82516SHong Zhang   MPI_Comm          comm;
9720d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
9738cb82516SHong Zhang #endif
9740d3441aeSHong Zhang 
9750d3441aeSHong Zhang   PetscFunctionBegin;
9760d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
9770d3441aeSHong Zhang 
978e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
9798cb82516SHong Zhang #if defined(PTAP_PROFILE)
9800d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
9818cb82516SHong Zhang #endif
982748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
9830d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
9840d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
985748c7196SHong Zhang   }
9868cb82516SHong Zhang #if defined(PTAP_PROFILE)
9870d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
9880d3441aeSHong Zhang   eR = t1 - t0;
9898cb82516SHong Zhang #endif
9900d3441aeSHong Zhang 
991e497d3c8SHong Zhang   /* 2) get AP_loc */
9920d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
9938cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
9940d3441aeSHong Zhang 
995e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
9960d3441aeSHong Zhang   /*-----------------------------------------------------*/
997748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
998748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
9990d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
10000d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1001e497d3c8SHong Zhang   }
10020d3441aeSHong Zhang 
10038cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
10048cb82516SHong Zhang   /* ---------------------------------------------- */
10050d3441aeSHong Zhang   /* get data from symbolic products */
10068cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
10072dd9e643SHong Zhang   if (ptap->P_oth) {
10088cb82516SHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
10092dd9e643SHong Zhang   }
10108cb82516SHong Zhang   apa   = ptap->apa;
1011681d504bSHong Zhang   api   = ap->i;
1012681d504bSHong Zhang   apj   = ap->j;
1013e497d3c8SHong Zhang   for (i=0; i<am; i++) {
1014445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
1015e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
1016e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
1017e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
1018e497d3c8SHong Zhang       col = apj[j+api[i]];
1019e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
1020e497d3c8SHong Zhang       apa[col] = 0.0;
10210d3441aeSHong Zhang     }
1022aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
10230d3441aeSHong Zhang   }
10248cb82516SHong Zhang #if defined(PTAP_PROFILE)
10250d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
10260d3441aeSHong Zhang   eAP = t2 - t1;
10278cb82516SHong Zhang #endif
10280d3441aeSHong Zhang 
10298cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1030445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1031445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
10329ce11a7cSHong Zhang   C_loc = ptap->C_loc;
10339ce11a7cSHong Zhang   C_oth = ptap->C_oth;
10348cb82516SHong Zhang #if defined(PTAP_PROFILE)
10350d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
10360d3441aeSHong Zhang   eCseq = t3 - t2;
10378cb82516SHong Zhang #endif
10380d3441aeSHong Zhang 
10390d3441aeSHong Zhang   /* add C_loc and Co to to C */
10400d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
10410d3441aeSHong Zhang 
10420d3441aeSHong Zhang   /* C_loc -> C */
10430d3441aeSHong Zhang   cm    = C_loc->rmap->N;
10449ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
10458cb82516SHong Zhang   cols = c_seq->j;
10468cb82516SHong Zhang   vals = c_seq->a;
10470d3441aeSHong Zhang   for (i=0; i<cm; i++) {
10489ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
10490d3441aeSHong Zhang     row = rstart + 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 
10540d3441aeSHong Zhang   /* Co -> C, off-processor part */
10559ce11a7cSHong Zhang   cm = C_oth->rmap->N;
10569ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
10578cb82516SHong Zhang   cols = c_seq->j;
10588cb82516SHong Zhang   vals = c_seq->a;
10599ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
10609ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
10610d3441aeSHong Zhang     row = p->garray[i];
10620d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10638cb82516SHong Zhang     cols += ncols; vals += ncols;
10640d3441aeSHong Zhang   }
10650d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10660d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10678cb82516SHong Zhang #if defined(PTAP_PROFILE)
10680d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
10690d3441aeSHong Zhang   eCmpi = t4 - t3;
10700d3441aeSHong Zhang 
10718cb82516SHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
10728cb82516SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
10730d3441aeSHong Zhang   if (rank==1) {
10740d3441aeSHong 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);
10750d3441aeSHong Zhang   }
10768cb82516SHong Zhang #endif
1077748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
10780d3441aeSHong Zhang   PetscFunctionReturn(0);
10790d3441aeSHong Zhang }
10800d3441aeSHong Zhang 
10818403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
108242c57489SHong Zhang {
108342c57489SHong Zhang   PetscErrorCode      ierr;
108477584fe6SHong Zhang   Mat                 Cmpi;
1085f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
10860298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1087f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
108842c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
108942c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1090ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
109177584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
1092a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1093d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
109442c57489SHong Zhang   PetscBT             lnkbt;
1095ce94432eSBarry Smith   MPI_Comm            comm;
1096a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
109742c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
109842c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
109924ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
110042c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1101ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
1102ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
110342c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
110477584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
1105d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
110642c57489SHong Zhang 
110742c57489SHong Zhang   PetscFunctionBegin;
1108ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
110983445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
111083445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
111183445d95SHong Zhang 
1112f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1113b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1114b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
11159f4155fbSHong Zhang   ptap->merge = merge;
1116f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
11176c6e00e0SHong Zhang 
11186c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1119b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1120fe615159SHong Zhang 
11216c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
1122a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
11236c6e00e0SHong Zhang 
1124a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1125a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
11266c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
11276c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
112842c57489SHong Zhang 
11292addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
11302addb229SHong Zhang   /*--------------------------------------------------------------------------*/
1131854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
113282412ba7SHong Zhang   api[0] = 0;
113383445d95SHong Zhang 
113483445d95SHong Zhang   /* create and initialize a linked list */
1135a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
113683445d95SHong Zhang 
1137a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
1138f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
1139f4a743e1SHong Zhang   current_space = free_space;
1140f4a743e1SHong Zhang 
1141f4a743e1SHong Zhang   for (i=0; i<am; i++) {
1142f4a743e1SHong Zhang     /* diagonal portion of A */
1143ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
114477584fe6SHong Zhang     aj  = ad->j + adi[i];
1145ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
114677584fe6SHong Zhang       row  = aj[j];
114782412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
1148ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
114983445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
1150a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1151f4a743e1SHong Zhang     }
1152f4a743e1SHong Zhang     /* off-diagonal portion of A */
1153ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
115477584fe6SHong Zhang     aj  = ao->j + aoi[i];
1155ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
115677584fe6SHong Zhang       row  = aj[j];
115782412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
1158ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
1159a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1160f4a743e1SHong Zhang     }
1161a914927fSHong Zhang     apnz     = lnk[0];
116282412ba7SHong Zhang     api[i+1] = api[i] + apnz;
116377584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
1164f4a743e1SHong Zhang 
116583445d95SHong Zhang     /* if free space is not available, double the total space in the list */
116682412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1167f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1168f2b054eeSHong Zhang       nspacedouble++;
1169f4a743e1SHong Zhang     }
1170f4a743e1SHong Zhang 
1171f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
1172a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11732205254eSKarl Rupp 
117482412ba7SHong Zhang     current_space->array           += apnz;
117582412ba7SHong Zhang     current_space->local_used      += apnz;
117682412ba7SHong Zhang     current_space->local_remaining -= apnz;
1177f4a743e1SHong Zhang   }
1178a914927fSHong Zhang 
117982412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
1180f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
1181854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
118277584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
1183118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
1184d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1185f4a743e1SHong Zhang 
11862addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
11872addb229SHong Zhang   /*-----------------------------------------------------------------*/
1188de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
118942c57489SHong Zhang 
1190ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
1191d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
119283445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1193854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1194de0260b3SHong Zhang   coi[0] = 0;
119542c57489SHong Zhang 
1196d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
1197f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am]));
119822d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
119942c57489SHong Zhang   current_space = free_space;
120042c57489SHong Zhang 
1201de0260b3SHong Zhang   for (i=0; i<pon; i++) {
120282412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
120377584fe6SHong Zhang     ptJ = potj + poti[i];
120477584fe6SHong Zhang     for (j=0; j<pnz; j++) {
120577584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
120682412ba7SHong Zhang       apnz = api[row+1] - api[row];
1207ded4f561SHong Zhang       Jptr = apj + api[row];
120883445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1209a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
121042c57489SHong Zhang     }
1211a914927fSHong Zhang     nnz = lnk[0];
121242c57489SHong Zhang 
121383445d95SHong Zhang     /* If free space is not available, double the total space in the list */
1214ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
1215f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1216d16ca5f0SHong Zhang       nspacedouble++;
121742c57489SHong Zhang     }
121842c57489SHong Zhang 
121942c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
1220a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
12212205254eSKarl Rupp 
1222ded4f561SHong Zhang     current_space->array           += nnz;
1223ded4f561SHong Zhang     current_space->local_used      += nnz;
1224ded4f561SHong Zhang     current_space->local_remaining -= nnz;
12252205254eSKarl Rupp 
1226ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
122742c57489SHong Zhang   }
12286b911d16SHong Zhang 
12296b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
1230a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1231118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
1232d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1233de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
123442c57489SHong Zhang 
12356b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
12366b911d16SHong Zhang   /*--------------------------------------------------*/
12376b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
12386b911d16SHong Zhang   len_s        = merge->len_s;
12396b911d16SHong Zhang   merge->nsend = 0;
12406b911d16SHong Zhang 
12416b911d16SHong Zhang 
124242c57489SHong Zhang   /* determine row ownership */
124326283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12447a2fc3feSBarry Smith   merge->rowmap->n  = pn;
12457a2fc3feSBarry Smith   merge->rowmap->bs = 1;
12462205254eSKarl Rupp 
124726283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
12487a2fc3feSBarry Smith   owners = merge->rowmap->range;
124942c57489SHong Zhang 
125042c57489SHong Zhang   /* determine the number of messages to send, their lengths */
1251dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
125283445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1253854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1254de0260b3SHong Zhang 
125583445d95SHong Zhang   proc = 0;
1256de0260b3SHong Zhang   for (i=0; i<pon; i++) {
1257de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
12582addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
1259ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
126083445d95SHong Zhang   }
1261de0260b3SHong Zhang 
12622addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
1263de0260b3SHong Zhang   owners_co[0] = 0;
126442c57489SHong Zhang   for (proc=0; proc<size; proc++) {
1265de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1266ee6e4b5bSHong Zhang     if (len_s[proc]) {
126742c57489SHong Zhang       merge->nsend++;
12682addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
126942c57489SHong Zhang       len         += len_si[proc];
127042c57489SHong Zhang     }
127142c57489SHong Zhang   }
127242c57489SHong Zhang 
1273ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12740298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
127542c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
127642c57489SHong Zhang 
1277ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
1278529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1279529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1280854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
128142c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
128242c57489SHong Zhang     if (!len_s[proc]) continue;
1283de0260b3SHong Zhang     i    = owners_co[proc];
1284529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
128542c57489SHong Zhang     k++;
128642c57489SHong Zhang   }
128742c57489SHong Zhang 
1288ded4f561SHong Zhang   /* receives and sends of coj are complete */
128977584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1290c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1291ded4f561SHong Zhang   }
1292ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
12930c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
129482412ba7SHong Zhang 
12956b911d16SHong Zhang   /* (4) send and recv coi */
12966b911d16SHong Zhang   /*-----------------------*/
1297529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1298529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1299854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
130042c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
130142c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
130242c57489SHong Zhang     if (!len_s[proc]) continue;
130342c57489SHong Zhang     /* form outgoing message for i-structure:
130442c57489SHong Zhang          buf_si[0]:                 nrows to be sent
130542c57489SHong Zhang                [1:nrows]:           row index (global)
130642c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
130742c57489SHong Zhang     */
130842c57489SHong Zhang     /*-------------------------------------------*/
13092addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
131042c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
131142c57489SHong Zhang     buf_si[0]   = nrows;
131242c57489SHong Zhang     buf_si_i[0] = 0;
131342c57489SHong Zhang     nrows       = 0;
1314de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1315ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
1316ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1317de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
131842c57489SHong Zhang       nrows++;
131942c57489SHong Zhang     }
1320529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
132142c57489SHong Zhang     k++;
132242c57489SHong Zhang     buf_si += len_si[proc];
132342c57489SHong Zhang   }
1324ded4f561SHong Zhang   i = merge->nrecv;
1325ded4f561SHong Zhang   while (i--) {
1326c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1327ded4f561SHong Zhang   }
1328ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
13290c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1330a914927fSHong Zhang 
133124ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
133242c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1333ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
133442c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
133542c57489SHong Zhang 
13366b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
13376b911d16SHong Zhang   /*----------------------------------------------*/
1338ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1339ded4f561SHong Zhang 
134024ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
1341854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
134224ecddacSHong Zhang   pti[0] = 0;
134342c57489SHong Zhang 
1344d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1345f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am]));
134622d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
134742c57489SHong Zhang   current_space = free_space;
134842c57489SHong Zhang 
1349dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
135042c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
135142c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
135242c57489SHong Zhang     nrows       = *buf_ri_k[k];
135342c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
135442c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
135542c57489SHong Zhang   }
135642c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
1357936ad1eaSHong Zhang 
135842c57489SHong Zhang   for (i=0; i<pn; i++) {
1359ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
1360ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
136177584fe6SHong Zhang     ptJ = pdtj + pdti[i];
136277584fe6SHong Zhang     for (j=0; j<pnz; j++) {
136377584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1364ded4f561SHong Zhang       apnz = api[row+1] - api[row];
1365ded4f561SHong Zhang       Jptr = apj + api[row];
1366ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1367a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1368ded4f561SHong Zhang     }
1369a914927fSHong Zhang 
137042c57489SHong Zhang     /* add received col data into lnk */
137142c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
137242c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1373ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1374ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1375a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
137642c57489SHong Zhang         nextrow[k]++; nextci[k]++;
137742c57489SHong Zhang       }
137842c57489SHong Zhang     }
1379a914927fSHong Zhang     nnz = lnk[0];
138042c57489SHong Zhang 
138142c57489SHong Zhang     /* if free space is not available, make more free space */
1382ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
1383f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1384d16ca5f0SHong Zhang       nspacedouble++;
138542c57489SHong Zhang     }
138642c57489SHong Zhang     /* copy data into free space, then initialize lnk */
1387a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1388ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13892205254eSKarl Rupp 
1390ded4f561SHong Zhang     current_space->array           += nnz;
1391ded4f561SHong Zhang     current_space->local_used      += nnz;
1392ded4f561SHong Zhang     current_space->local_remaining -= nnz;
13932205254eSKarl Rupp 
139424ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
139542c57489SHong Zhang   }
1396ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
13970572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
139842c57489SHong Zhang 
1399854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
140024ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
140124ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
1402d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
140342c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
140442c57489SHong Zhang 
14056b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
14066b911d16SHong Zhang   /*------------------------------------------*/
140777584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
140877584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
140933d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
141077584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
141177584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
141242c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1413a2f3521dSMark F. Adams 
1414b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
1415b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
1416b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
1417b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
141842c57489SHong Zhang   merge->buf_ri    = buf_ri;
141942c57489SHong Zhang   merge->buf_rj    = buf_rj;
1420de0260b3SHong Zhang   merge->owners_co = owners_co;
142177584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
142242c57489SHong Zhang 
142377584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
142477584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1425f8487c73SHong Zhang   c->ptap     = ptap;
142677584fe6SHong Zhang   ptap->api   = api;
142777584fe6SHong Zhang   ptap->apj   = apj;
14288403a639SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
14298403a639SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
14308403a639SHong Zhang 
14318403a639SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
14328403a639SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
14338403a639SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
14348403a639SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
14358403a639SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
143677584fe6SHong Zhang   *C                     = Cmpi;
1437414894bdSHong Zhang 
1438414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
1439414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
1440414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
1441414894bdSHong Zhang   /* set default scalable */
1442da0a95b2SSatish Balay   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
14432205254eSKarl Rupp 
1444c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
1445414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
14461795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1447414894bdSHong Zhang   } else {
14481795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
1449414894bdSHong Zhang   }
1450414894bdSHong Zhang 
1451f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
145224ecddacSHong Zhang   if (pti[pn] != 0) {
145357622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
145457622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1455f2b054eeSHong Zhang   } else {
145677584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1457f2b054eeSHong Zhang   }
1458f2b054eeSHong Zhang #endif
145942c57489SHong Zhang   PetscFunctionReturn(0);
146042c57489SHong Zhang }
146142c57489SHong Zhang 
14628403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
146342c57489SHong Zhang {
146442c57489SHong Zhang   PetscErrorCode      ierr;
1465f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
146642c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1467de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
146842c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1469f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
14709f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
14711d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
147282412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
147382412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1474e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1475d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1476ce94432eSBarry Smith   MPI_Comm            comm;
147742c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
147882412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
147942c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
148050f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
148142c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
148242c57489SHong Zhang   MPI_Status          *status;
148382412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
148482412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
1485d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
14866ce94e8fSHong Zhang   PetscBool           scalable;
148738571addSHong Zhang #if defined(PTAP_PROFILE)
1488e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
148938571addSHong Zhang #endif
149042c57489SHong Zhang 
149142c57489SHong Zhang   PetscFunctionBegin;
1492ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
149338571addSHong Zhang #if defined(PTAP_PROFILE)
14948563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
149538571addSHong Zhang #endif
149642c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
149742c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
149842c57489SHong Zhang 
1499f8487c73SHong Zhang   ptap = c->ptap;
1500ce94432eSBarry 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");
1501f8487c73SHong Zhang   merge    = ptap->merge;
1502414894bdSHong Zhang   apa      = ptap->apa;
15036ce94e8fSHong Zhang   scalable = ptap->scalable;
15046c6e00e0SHong Zhang 
1505a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
15066b911d16SHong Zhang   /*-----------------------------------------------------*/
1507f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
15089f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1509f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
15106c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1511b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1512a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
15136c6e00e0SHong Zhang   }
151438571addSHong Zhang #if defined(PTAP_PROFILE)
15158563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
151605d62848SHong Zhang   eP = t1-t0;
151738571addSHong Zhang #endif
151805d62848SHong Zhang   /*
151905d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
152005d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
152105d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
152205d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
152305d62848SHong Zhang    */
1524f8487c73SHong Zhang 
1525f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1526f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
152742c57489SHong Zhang   /* get data from symbolic products */
1528a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1529a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
153089ae1891SBarry Smith   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
153142c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
153242c57489SHong Zhang 
1533de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
15341795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1535de0260b3SHong Zhang 
1536de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
15377a2fc3feSBarry Smith   owners = merge->rowmap->range;
15381795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
153942c57489SHong Zhang 
1540a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1541d9824c15SHong Zhang 
15429516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1543b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
15448cb82516SHong Zhang #if defined(PTAP_PROFILE)
154505d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
15468cb82516SHong Zhang #endif
1547a9d06231SHong Zhang     for (i=0; i<am; i++) {
1548b091e509SHong Zhang #if defined(PTAP_PROFILE)
15498563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1550b091e509SHong Zhang #endif
1551a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1552a9d06231SHong Zhang       /*------------------------------------------------------------*/
1553a9d06231SHong Zhang       apJ = apj + api[i];
1554a9d06231SHong Zhang 
1555a9d06231SHong Zhang       /* diagonal portion of A */
1556a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1557a9d06231SHong Zhang       adj = ad->j + adi[i];
1558a9d06231SHong Zhang       ada = ad->a + adi[i];
1559a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1560a9d06231SHong Zhang         row = adj[j];
1561a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1562a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1563a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1564a9d06231SHong Zhang 
1565a9d06231SHong Zhang         /* perform dense axpy */
1566e900a757SHong Zhang         valtmp = ada[j];
1567a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1568e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1569a9d06231SHong Zhang         }
1570a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1571a9d06231SHong Zhang       }
1572a9d06231SHong Zhang 
1573a9d06231SHong Zhang       /* off-diagonal portion of A */
1574a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1575a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1576a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1577a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1578a9d06231SHong Zhang         row = aoj[j];
1579a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1580a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1581a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1582a9d06231SHong Zhang 
1583a9d06231SHong Zhang         /* perform dense axpy */
1584e900a757SHong Zhang         valtmp = aoa[j];
1585a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1586e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1587a9d06231SHong Zhang         }
1588a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1589a9d06231SHong Zhang       }
1590b091e509SHong Zhang #if defined(PTAP_PROFILE)
15918563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1592b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1593b091e509SHong Zhang #endif
1594a9d06231SHong Zhang 
1595a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1596a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1597a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1598a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1599a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1600e900a757SHong Zhang       poJ = po->j + po->i[i];
1601a9d06231SHong Zhang       pA  = po->a + po->i[i];
1602a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1603e900a757SHong Zhang         row = poJ[j];
1604e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1605e900a757SHong Zhang         cj  = coj + coi[row];
1606e900a757SHong Zhang         ca  = coa + coi[row];
1607a9d06231SHong Zhang         /* perform dense axpy */
1608e900a757SHong Zhang         valtmp = pA[j];
1609a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1610e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1611a9d06231SHong Zhang         }
1612a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1613a9d06231SHong Zhang       }
1614a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1615a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1616e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1617a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1618a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1619e900a757SHong Zhang         row = pdJ[j];
1620e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1621e900a757SHong Zhang         cj  = bj + bi[row];
1622e900a757SHong Zhang         ca  = ba + bi[row];
1623a9d06231SHong Zhang         /* perform dense axpy */
1624e900a757SHong Zhang         valtmp = pA[j];
1625a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1626e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1627a9d06231SHong Zhang         }
1628a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1629a9d06231SHong Zhang       }
16308403a639SHong Zhang 
1631a9d06231SHong Zhang       /* zero the current row of A*P */
1632a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1633b091e509SHong Zhang #if defined(PTAP_PROFILE)
16348563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
163505d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1636b091e509SHong Zhang #endif
1637a9d06231SHong Zhang     }
163838571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1639b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
164038571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1641a9d06231SHong Zhang     pA=pa_loc;
164242c57489SHong Zhang     for (i=0; i<am; i++) {
1643b091e509SHong Zhang #if defined(PTAP_PROFILE)
16448563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1645b091e509SHong Zhang #endif
1646f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
164782412ba7SHong Zhang       apnz = api[i+1] - api[i];
164882412ba7SHong Zhang       apJ  = apj + api[i];
164942c57489SHong Zhang       /* diagonal portion of A */
165082412ba7SHong Zhang       anz = adi[i+1] - adi[i];
16511d633602SHong Zhang       adj = ad->j + adi[i];
16521d633602SHong Zhang       ada = ad->a + adi[i];
165382412ba7SHong Zhang       for (j=0; j<anz; j++) {
16541d633602SHong Zhang         row    = adj[j];
165582412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
165682412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
165782412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1658e900a757SHong Zhang         valtmp = ada[j];
1659f4a743e1SHong Zhang         nextp  = 0;
166082412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
166182412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1662e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
166342c57489SHong Zhang           }
166442c57489SHong Zhang         }
1665dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
166642c57489SHong Zhang       }
166742c57489SHong Zhang       /* off-diagonal portion of A */
166882412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
16691d633602SHong Zhang       aoj = ao->j + aoi[i];
16701d633602SHong Zhang       aoa = ao->a + aoi[i];
167182412ba7SHong Zhang       for (j=0; j<anz; j++) {
16721d633602SHong Zhang         row    = aoj[j];
167382412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
167482412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
167582412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1676e900a757SHong Zhang         valtmp = aoa[j];
1677f4a743e1SHong Zhang         nextp  = 0;
167882412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
167982412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1680e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
168142c57489SHong Zhang           }
168242c57489SHong Zhang         }
1683dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
168442c57489SHong Zhang       }
1685b091e509SHong Zhang #if defined(PTAP_PROFILE)
16868563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1687b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1688b091e509SHong Zhang #endif
168942c57489SHong Zhang 
1690a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1691a9d06231SHong Zhang       /*--------------------------------------------------------------*/
169282412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1693f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
169482412ba7SHong Zhang       for (j=0; j<pnz; j++) {
169542c57489SHong Zhang         nextap = 0;
1696f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
169782412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1698e900a757SHong Zhang           row = *poJ;
1699e900a757SHong Zhang           cj  = coj + coi[row];
1700e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
170182412ba7SHong Zhang         } else {                            /* put the value into Cd */
1702e900a757SHong Zhang           row = *pdJ;
1703e900a757SHong Zhang           cj  = bj + bi[row];
1704e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
170542c57489SHong Zhang         }
1706e900a757SHong Zhang         valtmp = pA[j];
170782412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1708e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
170942c57489SHong Zhang         }
1710dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1711de0260b3SHong Zhang       }
17120496f32cSHong Zhang       pA += pnz;
171342c57489SHong Zhang       /* zero the current row info for A*P */
171482412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1715b091e509SHong Zhang #if defined(PTAP_PROFILE)
17168563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
171705d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1718b091e509SHong Zhang #endif
171942c57489SHong Zhang     }
1720d9824c15SHong Zhang   }
172138571addSHong Zhang #if defined(PTAP_PROFILE)
17228563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
172338571addSHong Zhang #endif
172442c57489SHong Zhang 
1725a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1726a9d06231SHong Zhang   /*------------------------------------*/
172742c57489SHong Zhang   buf_ri = merge->buf_ri;
172842c57489SHong Zhang   buf_rj = merge->buf_rj;
172942c57489SHong Zhang   len_s  = merge->len_s;
1730fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
173142c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
173242c57489SHong Zhang 
1733dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
173442c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
173542c57489SHong Zhang     if (!len_s[proc]) continue;
1736de0260b3SHong Zhang     i    = merge->owners_co[proc];
1737de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
173842c57489SHong Zhang     k++;
173942c57489SHong Zhang   }
17400c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
17410c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
174242c57489SHong Zhang 
1743a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
174442c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
174582412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
174638571addSHong Zhang #if defined(PTAP_PROFILE)
17478563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
174838571addSHong Zhang #endif
174942c57489SHong Zhang 
1750a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1751a9d06231SHong Zhang   /*------------------------------------------------------*/
1752dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
175342c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
175442c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
175542c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
175642c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
175782412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
175842c57489SHong Zhang   }
175942c57489SHong Zhang 
1760de0260b3SHong Zhang   for (i=0; i<cm; i++) {
176182412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
176242c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1763de0260b3SHong Zhang     ba_i = ba + bi[i];
176482412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
176542c57489SHong Zhang     /* add received vals into ba */
176642c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
176742c57489SHong Zhang       /* i-th row */
176842c57489SHong Zhang       if (i == *nextrow[k]) {
176982412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
177082412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
177182412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
177282412ba7SHong Zhang         nextcj = 0;
177382412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
177482412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
177582412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
177642c57489SHong Zhang           }
177742c57489SHong Zhang         }
177882412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1779c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
178042c57489SHong Zhang       }
178142c57489SHong Zhang     }
178282412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
178342c57489SHong Zhang   }
178442c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
178542c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
178642c57489SHong Zhang 
1787de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1788c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
178942c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
17900572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
179138571addSHong Zhang #if defined(PTAP_PROFILE)
17928563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
17939516a85cSHong Zhang   if (rank==1) {
179405d62848SHong 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);
17959516a85cSHong Zhang   }
179638571addSHong Zhang #endif
179742c57489SHong Zhang   PetscFunctionReturn(0);
179842c57489SHong Zhang }
1799