xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 904d1e7050d8278cb5f67df70e5edfa4eaab048b)
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 
15cf97cf3cSHong Zhang PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer)
16cf97cf3cSHong Zhang {
17cf97cf3cSHong Zhang   PetscErrorCode    ierr;
18cf97cf3cSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data;
19cf97cf3cSHong Zhang   Mat_PtAPMPI       *ptap=a->ptap;
20cf97cf3cSHong Zhang   PetscBool         iascii;
21cf97cf3cSHong Zhang   PetscViewerFormat format;
22cf97cf3cSHong Zhang 
23cf97cf3cSHong Zhang   PetscFunctionBegin;
24cf97cf3cSHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
25cf97cf3cSHong Zhang   if (iascii) {
26cf97cf3cSHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
27cf97cf3cSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
28cf97cf3cSHong Zhang       if (ptap->algType == 0) {
29cf97cf3cSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr);
30cf97cf3cSHong Zhang       } else if (ptap->algType == 1) {
31cf97cf3cSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr);
32cf97cf3cSHong Zhang       }
33cf97cf3cSHong Zhang     }
34cf97cf3cSHong Zhang   }
35cf97cf3cSHong Zhang   ierr = (ptap->view)(A,viewer);CHKERRQ(ierr);
36cf97cf3cSHong Zhang   PetscFunctionReturn(0);
37cf97cf3cSHong Zhang }
38cf97cf3cSHong Zhang 
39f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
40cc6cb787SHong Zhang {
41cc6cb787SHong Zhang   PetscErrorCode ierr;
42f8487c73SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
43f8487c73SHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
44cc6cb787SHong Zhang 
45cc6cb787SHong Zhang   PetscFunctionBegin;
46f8487c73SHong Zhang   if (ptap) {
47c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI *merge=ptap->merge;
48b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
49f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
50a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
51a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
52c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
5305d62848SHong Zhang     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
5405d62848SHong Zhang     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
558403a639SHong Zhang     if (ptap->AP_loc) { /* used by alg_rap */
56681d504bSHong Zhang       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
57681d504bSHong Zhang       ierr = PetscFree(ap->i);CHKERRQ(ierr);
58681d504bSHong Zhang       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
590d3441aeSHong Zhang       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
608403a639SHong Zhang     } else { /* used by alg_ptap */
618403a639SHong Zhang       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
628403a639SHong Zhang       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
63681d504bSHong Zhang     }
642259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
652259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
66414894bdSHong Zhang     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
67a560ef98SHong Zhang 
688403a639SHong Zhang     if (merge) { /* used by alg_ptap */
69cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
70cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
71cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
72cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
73cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
74c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
75cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
76c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
77cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
78445158ffSHong Zhang       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
79445158ffSHong Zhang       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
8005b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
816bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
82f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
83bf0cc555SLisandro Dalcin     }
84a560ef98SHong Zhang 
85a560ef98SHong Zhang     ierr = ptap->destroy(A);CHKERRQ(ierr);
86f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
87c8b0795cSMark F. Adams   }
88cc6cb787SHong Zhang   PetscFunctionReturn(0);
89cc6cb787SHong Zhang }
90cc6cb787SHong Zhang 
91aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
921a47ec54SHong Zhang {
931a47ec54SHong Zhang   PetscErrorCode ierr;
941a47ec54SHong Zhang   Mat_MPIAIJ     *a     = (Mat_MPIAIJ*)A->data;
951a47ec54SHong Zhang   Mat_PtAPMPI    *ptap  = a->ptap;
961a47ec54SHong Zhang 
971a47ec54SHong Zhang   PetscFunctionBegin;
981a47ec54SHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
991a47ec54SHong Zhang   (*M)->ops->destroy   = ptap->destroy;
1001a47ec54SHong Zhang   (*M)->ops->duplicate = ptap->duplicate;
101cff31925SHong Zhang   (*M)->ops->view      = ptap->view;
1021a47ec54SHong Zhang   PetscFunctionReturn(0);
1031a47ec54SHong Zhang }
1041a47ec54SHong Zhang 
105150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
10642c57489SHong Zhang {
10742c57489SHong Zhang   PetscErrorCode ierr;
108a4ffb656SHong Zhang   PetscBool      flg;
10967a12041SHong Zhang   MPI_Comm       comm;
110a4ffb656SHong Zhang #if !defined(PETSC_HAVE_HYPRE)
111a4ffb656SHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
112a4ffb656SHong Zhang   PetscInt            nalg=2;
113a4ffb656SHong Zhang #else
114a4ffb656SHong Zhang   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
115a4ffb656SHong Zhang   PetscInt            nalg=3;
116a4ffb656SHong Zhang #endif
117a4ffb656SHong Zhang   PetscInt            pN=P->cmap->N,alg=1; /* set default algorithm */
11842c57489SHong Zhang 
11942c57489SHong Zhang   PetscFunctionBegin;
12067a12041SHong Zhang   /* check if matrix local sizes are compatible */
12167a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1226c4ed002SBarry 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);
1236c4ed002SBarry 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);
12467a12041SHong Zhang 
125cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
126a4ffb656SHong Zhang     /* pick an algorithm */
127a4ffb656SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
128a4ffb656SHong Zhang     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
129a4ffb656SHong Zhang     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
130a4ffb656SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
131a4ffb656SHong Zhang 
132a4ffb656SHong Zhang     if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
133a4ffb656SHong Zhang       MatInfo     Ainfo,Pinfo;
134a4ffb656SHong Zhang       PetscInt    nz_local;
135a4ffb656SHong Zhang       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
136a4ffb656SHong Zhang 
137a4ffb656SHong Zhang       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
138a4ffb656SHong Zhang       ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
139a4ffb656SHong Zhang       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
140a4ffb656SHong Zhang 
141a4ffb656SHong Zhang       if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
142a4ffb656SHong Zhang       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
143a4ffb656SHong Zhang 
144a4ffb656SHong Zhang       if (alg_scalable) {
145a4ffb656SHong Zhang         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
1460d3441aeSHong Zhang       }
147a4ffb656SHong Zhang     }
148a4ffb656SHong Zhang 
149a4ffb656SHong Zhang     switch (alg) {
150a4ffb656SHong Zhang     case 1:
151a4ffb656SHong Zhang       /* do R=P^T locally, then C=R*A*P */
152a4ffb656SHong Zhang       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
153a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
1543ff4c91cSHong Zhang       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
155a4ffb656SHong Zhang       break;
156a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE)
157a4ffb656SHong Zhang     case 2:
158a4ffb656SHong Zhang       /* Use boomerAMGBuildCoarseOperator */
159a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
160a4ffb656SHong Zhang       PetscFunctionReturn(0);
161a4ffb656SHong Zhang       break;
162a4ffb656SHong Zhang #endif
163a4ffb656SHong Zhang     default:
164a4ffb656SHong Zhang       /* do R=P^T locally, then C=R*A*P */
165a4ffb656SHong Zhang       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
166a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
167a4ffb656SHong Zhang       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
168a4ffb656SHong Zhang       break;
169a4ffb656SHong Zhang     }
1707a7894deSKris Buschelman   }
1713ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1728403a639SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
1733ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
17442c57489SHong Zhang   PetscFunctionReturn(0);
17542c57489SHong Zhang }
17642c57489SHong Zhang 
177cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
178cf742fccSHong Zhang {
179cf742fccSHong Zhang   PetscErrorCode    ierr;
180cf742fccSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
181cf742fccSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
182cf742fccSHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
183cf742fccSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
184cf742fccSHong Zhang   Mat               AP_loc,C_loc,C_oth;
1855ca39647SHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz;
186cf742fccSHong Zhang   PetscScalar       *apa;
187cf742fccSHong Zhang   const PetscInt    *cols;
188cf742fccSHong Zhang   const PetscScalar *vals;
189cf742fccSHong Zhang 
190cf742fccSHong Zhang   PetscFunctionBegin;
191cf742fccSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
192cf742fccSHong Zhang 
193cf742fccSHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
194cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
195cf742fccSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
196cf742fccSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
197cf742fccSHong Zhang   }
198cf742fccSHong Zhang 
199cf742fccSHong Zhang   /* 2) get AP_loc */
200cf742fccSHong Zhang   AP_loc = ptap->AP_loc;
201cf742fccSHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
202cf742fccSHong Zhang 
203cf742fccSHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
204cf742fccSHong Zhang   /*-----------------------------------------------------*/
205cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
206cf742fccSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
207cf742fccSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
208cf742fccSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
209cf742fccSHong Zhang   }
210cf742fccSHong Zhang 
211cf742fccSHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
212cf742fccSHong Zhang   /* ---------------------------------------------- */
213cf742fccSHong Zhang   /* get data from symbolic products */
214cf742fccSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
21552f7967eSHong Zhang   if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
216c072d3e2SSatish Balay   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!");
21752f7967eSHong Zhang 
218cf742fccSHong Zhang   api   = ap->i;
219cf742fccSHong Zhang   apj   = ap->j;
220cf742fccSHong Zhang   for (i=0; i<am; i++) {
221cf742fccSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
222cf742fccSHong Zhang     apnz = api[i+1] - api[i];
223b4e8d1b6SHong Zhang     apa = ap->a + api[i];
224b4e8d1b6SHong Zhang     ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
225b4e8d1b6SHong Zhang     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
226cf742fccSHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
227cf742fccSHong Zhang   }
228cf742fccSHong Zhang 
229cf742fccSHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
230cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
231cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
232cf742fccSHong Zhang 
233cf742fccSHong Zhang   C_loc = ptap->C_loc;
234cf742fccSHong Zhang   C_oth = ptap->C_oth;
235cf742fccSHong Zhang 
236cf742fccSHong Zhang   /* add C_loc and Co to to C */
237cf742fccSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
238cf742fccSHong Zhang 
239cf742fccSHong Zhang   /* C_loc -> C */
240cf742fccSHong Zhang   cm    = C_loc->rmap->N;
241cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
242cf742fccSHong Zhang   cols = c_seq->j;
243cf742fccSHong Zhang   vals = c_seq->a;
244cf742fccSHong Zhang   for (i=0; i<cm; i++) {
245cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
246cf742fccSHong Zhang     row = rstart + i;
247cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
248cf742fccSHong Zhang     cols += ncols; vals += ncols;
249cf742fccSHong Zhang   }
250cf742fccSHong Zhang 
251cf742fccSHong Zhang   /* Co -> C, off-processor part */
252cf742fccSHong Zhang   cm = C_oth->rmap->N;
253cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
254cf742fccSHong Zhang   cols = c_seq->j;
255cf742fccSHong Zhang   vals = c_seq->a;
256cf742fccSHong Zhang   for (i=0; i<cm; i++) {
257cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
258cf742fccSHong Zhang     row = p->garray[i];
259cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
260cf742fccSHong Zhang     cols += ncols; vals += ncols;
261cf742fccSHong Zhang   }
262cf742fccSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
263cf742fccSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264cf742fccSHong Zhang 
265cf742fccSHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
266cf742fccSHong Zhang   PetscFunctionReturn(0);
267cf742fccSHong Zhang }
268cf742fccSHong Zhang 
269e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
2700d3441aeSHong Zhang {
2710d3441aeSHong Zhang   PetscErrorCode      ierr;
2720d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
273681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
2742259aa2eSHong Zhang   MPI_Comm            comm;
2752259aa2eSHong Zhang   PetscMPIInt         size,rank;
27676db11f6SHong Zhang   Mat                 Cmpi,P_loc,P_oth;
27715a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
278d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
279d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
280f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
28115a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
282d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
28315a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
28415a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
28515a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
286445158ffSHong Zhang   PetscLayout         rowmap;
287445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
288445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
289b4e8d1b6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi;
29052f7967eSHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
29167a12041SHong Zhang   PetscScalar         *apv;
29278d30b94SHong Zhang   PetscTable          ta;
293aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
2948cb82516SHong Zhang   PetscReal           apfill;
295aa690a28SHong Zhang #endif
29667a12041SHong Zhang 
29767a12041SHong Zhang   PetscFunctionBegin;
29867a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
29967a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
30067a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
301ae5f0867Sstefano_zampini 
30252f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
30352f7967eSHong Zhang 
304ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
305ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
306ae5f0867Sstefano_zampini   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
307ae5f0867Sstefano_zampini 
308e953e47cSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
309e953e47cSHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
310e953e47cSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
311cf97cf3cSHong Zhang   ptap->algType = 0;
312e953e47cSHong Zhang 
313e953e47cSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
31476db11f6SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
315e953e47cSHong Zhang   /* get P_loc by taking all local rows of P */
31676db11f6SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
31776db11f6SHong Zhang 
31876db11f6SHong Zhang   ptap->P_loc = P_loc;
31976db11f6SHong Zhang   ptap->P_oth = P_oth;
320e953e47cSHong Zhang 
321e953e47cSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
322e953e47cSHong Zhang   /* --------------------------------- */
323e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
324e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
325e953e47cSHong Zhang 
326e953e47cSHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
327e953e47cSHong Zhang   /* ----------------------------------------------------------------- */
32876db11f6SHong Zhang   p_loc  = (Mat_SeqAIJ*)P_loc->data;
32952f7967eSHong Zhang   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
330e953e47cSHong Zhang 
331e953e47cSHong Zhang   /* create and initialize a linked list */
332e953e47cSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
33376db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
33476db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
335e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
336e953e47cSHong Zhang 
33776db11f6SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
338e953e47cSHong Zhang 
339e953e47cSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
34052f7967eSHong Zhang   if (ao) {
341e953e47cSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
34252f7967eSHong Zhang   } else {
34352f7967eSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
34452f7967eSHong Zhang   }
345e953e47cSHong Zhang   current_space = free_space;
346e953e47cSHong Zhang   nspacedouble  = 0;
347e953e47cSHong Zhang 
348e953e47cSHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
349e953e47cSHong Zhang   api[0] = 0;
350e953e47cSHong Zhang   for (i=0; i<am; i++) {
351e953e47cSHong Zhang     /* diagonal portion: Ad[i,:]*P */
352e953e47cSHong Zhang     ai = ad->i; pi = p_loc->i;
353e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
354e953e47cSHong Zhang     aj  = ad->j + ai[i];
355e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
356e953e47cSHong Zhang       row  = aj[j];
357e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
358e953e47cSHong Zhang       Jptr = p_loc->j + pi[row];
359e953e47cSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
36076db11f6SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
361e953e47cSHong Zhang     }
362e953e47cSHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
36352f7967eSHong Zhang     if (ao) {
364e953e47cSHong Zhang       ai = ao->i; pi = p_oth->i;
365e953e47cSHong Zhang       nzi = ai[i+1] - ai[i];
366e953e47cSHong Zhang       aj  = ao->j + ai[i];
367e953e47cSHong Zhang       for (j=0; j<nzi; j++) {
368e953e47cSHong Zhang         row  = aj[j];
369e953e47cSHong Zhang         pnz  = pi[row+1] - pi[row];
370e953e47cSHong Zhang         Jptr = p_oth->j + pi[row];
37176db11f6SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
372e953e47cSHong Zhang       }
37352f7967eSHong Zhang     }
374e953e47cSHong Zhang     apnz     = lnk[0];
375e953e47cSHong Zhang     api[i+1] = api[i] + apnz;
376e953e47cSHong Zhang 
377e953e47cSHong Zhang     /* if free space is not available, double the total space in the list */
378e953e47cSHong Zhang     if (current_space->local_remaining<apnz) {
379e953e47cSHong Zhang       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
380e953e47cSHong Zhang       nspacedouble++;
381e953e47cSHong Zhang     }
382e953e47cSHong Zhang 
383e953e47cSHong Zhang     /* Copy data into free space, then initialize lnk */
38476db11f6SHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
385e953e47cSHong Zhang 
386e953e47cSHong Zhang     current_space->array           += apnz;
387e953e47cSHong Zhang     current_space->local_used      += apnz;
388e953e47cSHong Zhang     current_space->local_remaining -= apnz;
389e953e47cSHong Zhang   }
390e953e47cSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
391e953e47cSHong Zhang   /* destroy list of free space and other temporary array(s) */
392e953e47cSHong Zhang   ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
393e953e47cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
39476db11f6SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
395e953e47cSHong Zhang 
396e953e47cSHong Zhang   /* Create AP_loc for reuse */
397e953e47cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
398e953e47cSHong Zhang 
399e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
40052f7967eSHong Zhang   if (ao) {
401e953e47cSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
40252f7967eSHong Zhang   } else {
40352f7967eSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
40452f7967eSHong Zhang   }
405e953e47cSHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
406e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
407e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
408e953e47cSHong Zhang 
409e953e47cSHong Zhang   if (api[am]) {
410b11c1ec8SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"Scalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
411e953e47cSHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
412e953e47cSHong Zhang   } else {
413b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
414e953e47cSHong Zhang   }
415e953e47cSHong Zhang #endif
416e953e47cSHong Zhang 
417e953e47cSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
418e953e47cSHong Zhang   /* ------------------------------------ */
419e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
420e953e47cSHong Zhang 
421e953e47cSHong Zhang   /* (3) send coj of C_oth to other processors  */
422e953e47cSHong Zhang   /* ------------------------------------------ */
423e953e47cSHong Zhang   /* determine row ownership */
424e953e47cSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
425e953e47cSHong Zhang   rowmap->n  = pn;
426e953e47cSHong Zhang   rowmap->bs = 1;
427e953e47cSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
428e953e47cSHong Zhang   owners = rowmap->range;
429e953e47cSHong Zhang 
430e953e47cSHong Zhang   /* determine the number of messages to send, their lengths */
431e953e47cSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
432e953e47cSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
433e953e47cSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
434e953e47cSHong Zhang 
435e953e47cSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
436e953e47cSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
437e953e47cSHong Zhang   con   = ptap->C_oth->rmap->n;
438e953e47cSHong Zhang   proc  = 0;
439e953e47cSHong Zhang   for (i=0; i<con; i++) {
440e953e47cSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
441e953e47cSHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
442e953e47cSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
443e953e47cSHong Zhang   }
444e953e47cSHong Zhang 
445e953e47cSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
446e953e47cSHong Zhang   owners_co[0] = 0;
447e953e47cSHong Zhang   nsend        = 0;
448e953e47cSHong Zhang   for (proc=0; proc<size; proc++) {
449e953e47cSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
450e953e47cSHong Zhang     if (len_s[proc]) {
451e953e47cSHong Zhang       nsend++;
452e953e47cSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
453e953e47cSHong Zhang       len         += len_si[proc];
454e953e47cSHong Zhang     }
455e953e47cSHong Zhang   }
456e953e47cSHong Zhang 
457e953e47cSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
458e953e47cSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
459e953e47cSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
460e953e47cSHong Zhang 
461e953e47cSHong Zhang   /* post the Irecv and Isend of coj */
462e953e47cSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
463e953e47cSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
464e953e47cSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
465e953e47cSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
466e953e47cSHong Zhang     if (!len_s[proc]) continue;
467e953e47cSHong Zhang     i    = owners_co[proc];
468e953e47cSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
469e953e47cSHong Zhang     k++;
470e953e47cSHong Zhang   }
471e953e47cSHong Zhang 
472e953e47cSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
473e953e47cSHong Zhang   /* ---------------------------------------- */
474e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
475e953e47cSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
476e953e47cSHong Zhang 
477e953e47cSHong Zhang   /* receives coj are complete */
478e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
479e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
480e953e47cSHong Zhang   }
481e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
482e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
483e953e47cSHong Zhang 
484e953e47cSHong Zhang   /* add received column indices into ta to update Crmax */
485e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
486e953e47cSHong Zhang     Jptr = buf_rj[k];
487e953e47cSHong Zhang     for (j=0; j<len_r[k]; j++) {
488e953e47cSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
489e953e47cSHong Zhang     }
490e953e47cSHong Zhang   }
491e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
492e953e47cSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
493e953e47cSHong Zhang 
494e953e47cSHong Zhang   /* (4) send and recv coi */
495e953e47cSHong Zhang   /*-----------------------*/
496e953e47cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
497e953e47cSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
498e953e47cSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
499e953e47cSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
500e953e47cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
501e953e47cSHong Zhang     if (!len_s[proc]) continue;
502e953e47cSHong Zhang     /* form outgoing message for i-structure:
503e953e47cSHong Zhang          buf_si[0]:                 nrows to be sent
504e953e47cSHong Zhang                [1:nrows]:           row index (global)
505e953e47cSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
506e953e47cSHong Zhang     */
507e953e47cSHong Zhang     /*-------------------------------------------*/
508e953e47cSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
509e953e47cSHong Zhang     buf_si_i    = buf_si + nrows+1;
510e953e47cSHong Zhang     buf_si[0]   = nrows;
511e953e47cSHong Zhang     buf_si_i[0] = 0;
512e953e47cSHong Zhang     nrows       = 0;
513e953e47cSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
514e953e47cSHong Zhang       nzi = coi[i+1] - coi[i];
515e953e47cSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
516e953e47cSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
517e953e47cSHong Zhang       nrows++;
518e953e47cSHong Zhang     }
519e953e47cSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
520e953e47cSHong Zhang     k++;
521e953e47cSHong Zhang     buf_si += len_si[proc];
522e953e47cSHong Zhang   }
523e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
524e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
525e953e47cSHong Zhang   }
526e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
527e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
528e953e47cSHong Zhang 
529e953e47cSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
530e953e47cSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
531e953e47cSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
532e953e47cSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
533b4e8d1b6SHong Zhang 
534e953e47cSHong Zhang   /* (5) compute the local portion of Cmpi      */
535e953e47cSHong Zhang   /* ------------------------------------------ */
536e953e47cSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
537e953e47cSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
538e953e47cSHong Zhang   current_space = free_space;
539e953e47cSHong Zhang 
540e953e47cSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
541e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {
542e953e47cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
543e953e47cSHong Zhang     nrows       = *buf_ri_k[k];
544e953e47cSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
545e953e47cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
546e953e47cSHong Zhang   }
547e953e47cSHong Zhang 
548e953e47cSHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
549cf742fccSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
550e953e47cSHong Zhang   for (i=0; i<pn; i++) {
551e953e47cSHong Zhang     /* add C_loc into Cmpi */
552e953e47cSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
553e953e47cSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
554cf742fccSHong Zhang     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
555e953e47cSHong Zhang 
556e953e47cSHong Zhang     /* add received col data into lnk */
557e953e47cSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
558e953e47cSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
559e953e47cSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
560e953e47cSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
561cf742fccSHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
562e953e47cSHong Zhang         nextrow[k]++; nextci[k]++;
563e953e47cSHong Zhang       }
564e953e47cSHong Zhang     }
565e953e47cSHong Zhang     nzi = lnk[0];
566e953e47cSHong Zhang 
567e953e47cSHong Zhang     /* copy data into free space, then initialize lnk */
568cf742fccSHong Zhang     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
569e953e47cSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
570e953e47cSHong Zhang   }
571e953e47cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
572cf742fccSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
573e953e47cSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
574e953e47cSHong Zhang 
575e953e47cSHong Zhang   /* local sizes and preallocation */
576e953e47cSHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
577e953e47cSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
578e953e47cSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
579e953e47cSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
580e953e47cSHong Zhang 
581e953e47cSHong Zhang   /* members in merge */
582e953e47cSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
583e953e47cSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
584e953e47cSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
585e953e47cSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
586e953e47cSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
587e953e47cSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
588e953e47cSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
589e953e47cSHong Zhang 
590e953e47cSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
591e953e47cSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
592e953e47cSHong Zhang   c->ptap         = ptap;
593e953e47cSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
594e953e47cSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
595cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
596e953e47cSHong Zhang 
597e953e47cSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
598e953e47cSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
599a4ffb656SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
600e953e47cSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
601e953e47cSHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
602cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
603e953e47cSHong Zhang   *C                     = Cmpi;
604e953e47cSHong Zhang   PetscFunctionReturn(0);
605e953e47cSHong Zhang }
606e953e47cSHong Zhang 
607e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
608e953e47cSHong Zhang {
609e953e47cSHong Zhang   PetscErrorCode      ierr;
610e953e47cSHong Zhang   Mat_PtAPMPI         *ptap;
611e953e47cSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
612e953e47cSHong Zhang   MPI_Comm            comm;
613e953e47cSHong Zhang   PetscMPIInt         size,rank;
614e953e47cSHong Zhang   Mat                 Cmpi;
615e953e47cSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
616e953e47cSHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
617e953e47cSHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
618e953e47cSHong Zhang   PetscBT             lnkbt;
619e953e47cSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
620e953e47cSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
621e953e47cSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
622e953e47cSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
623e953e47cSHong Zhang   MPI_Request         *swaits,*rwaits;
624e953e47cSHong Zhang   MPI_Status          *sstatus,rstatus;
625e953e47cSHong Zhang   PetscLayout         rowmap;
626e953e47cSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
627e953e47cSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
628e953e47cSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
629ec07b8f8SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
630e953e47cSHong Zhang   PetscScalar         *apv;
631e953e47cSHong Zhang   PetscTable          ta;
632e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
633e953e47cSHong Zhang   PetscReal           apfill;
634e953e47cSHong Zhang #endif
635e953e47cSHong Zhang 
636e953e47cSHong Zhang   PetscFunctionBegin;
637b4e8d1b6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
638e953e47cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
639e953e47cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
640e953e47cSHong Zhang 
64152f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
642ec07b8f8SHong Zhang 
643e953e47cSHong Zhang   /* create symbolic parallel matrix Cmpi */
644e953e47cSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
645e953e47cSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
646e953e47cSHong Zhang 
647e953e47cSHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
648e953e47cSHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
649e953e47cSHong Zhang 
65015a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
65115a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
65215a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
653a4ffb656SHong Zhang   ptap->algType = 1;
65415a3b8e2SHong Zhang 
65515a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
65615a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
65715a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
65815a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
65915a3b8e2SHong Zhang 
66067a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
66167a12041SHong Zhang   /* --------------------------------- */
662de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
663de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
66415a3b8e2SHong Zhang 
66567a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
66667a12041SHong Zhang   /* ----------------------------------------------------------------- */
66767a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
66852f7967eSHong Zhang   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
669445158ffSHong Zhang 
6709a6dcab7SHong Zhang   /* create and initialize a linked list */
67145d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
6724b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
6734b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
67478d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
675d0e9a2f7SHong 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); */
67678d30b94SHong Zhang 
67778d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
678445158ffSHong Zhang 
6798cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
680ec07b8f8SHong Zhang   if (ao) {
681f91af8c7SBarry Smith     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
682ec07b8f8SHong Zhang   } else {
683ec07b8f8SHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
684ec07b8f8SHong Zhang   }
685445158ffSHong Zhang   current_space = free_space;
68667a12041SHong Zhang   nspacedouble  = 0;
68767a12041SHong Zhang 
68867a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
68967a12041SHong Zhang   api[0] = 0;
690445158ffSHong Zhang   for (i=0; i<am; i++) {
69167a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
69267a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
69367a12041SHong Zhang     nzi = ai[i+1] - ai[i];
69467a12041SHong Zhang     aj  = ad->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_loc->j + pi[row];
699445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
700445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
701445158ffSHong Zhang     }
70267a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
703ec07b8f8SHong Zhang     if (ao) {
70467a12041SHong Zhang       ai = ao->i; pi = p_oth->i;
70567a12041SHong Zhang       nzi = ai[i+1] - ai[i];
70667a12041SHong Zhang       aj  = ao->j + ai[i];
707445158ffSHong Zhang       for (j=0; j<nzi; j++) {
708445158ffSHong Zhang         row  = aj[j];
70967a12041SHong Zhang         pnz  = pi[row+1] - pi[row];
71067a12041SHong Zhang         Jptr = p_oth->j + pi[row];
711445158ffSHong Zhang         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
712445158ffSHong Zhang       }
713ec07b8f8SHong Zhang     }
714445158ffSHong Zhang     apnz     = lnk[0];
715445158ffSHong Zhang     api[i+1] = api[i] + apnz;
716445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
717445158ffSHong Zhang 
718445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
719445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
720f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
721445158ffSHong Zhang       nspacedouble++;
722445158ffSHong Zhang     }
723445158ffSHong Zhang 
724445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
725445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
726445158ffSHong Zhang 
727445158ffSHong Zhang     current_space->array           += apnz;
728445158ffSHong Zhang     current_space->local_used      += apnz;
729445158ffSHong Zhang     current_space->local_remaining -= apnz;
730445158ffSHong Zhang   }
731681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
732445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
733445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
734445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
7359a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
7369a6dcab7SHong Zhang 
737aa690a28SHong Zhang   /* Create AP_loc for reuse */
738445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
739aa690a28SHong Zhang 
740aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
741ec07b8f8SHong Zhang   if (ao) {
742aa690a28SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
743ec07b8f8SHong Zhang   } else {
744ec07b8f8SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
745ec07b8f8SHong Zhang   }
746aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
747aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
748aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
749aa690a28SHong Zhang 
750aa690a28SHong Zhang   if (api[am]) {
751b11c1ec8SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"Nonscalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
752aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
753aa690a28SHong Zhang   } else {
754b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
755aa690a28SHong Zhang   }
756aa690a28SHong Zhang #endif
757aa690a28SHong Zhang 
758681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
759681d504bSHong Zhang   /* ------------------------------------ */
76067a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
76115a3b8e2SHong Zhang 
76267a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
76367a12041SHong Zhang   /* ------------------------------------------ */
76415a3b8e2SHong Zhang   /* determine row ownership */
765445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
766445158ffSHong Zhang   rowmap->n  = pn;
767445158ffSHong Zhang   rowmap->bs = 1;
768445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
769445158ffSHong Zhang   owners = rowmap->range;
77015a3b8e2SHong Zhang 
77115a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
7728cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
7738cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
77415a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
77515a3b8e2SHong Zhang 
77667a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
77767a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
77867a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
77915a3b8e2SHong Zhang   proc  = 0;
78067a12041SHong Zhang   for (i=0; i<con; i++) {
78115a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
78215a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
78315a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
78415a3b8e2SHong Zhang   }
78515a3b8e2SHong Zhang 
78615a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
78715a3b8e2SHong Zhang   owners_co[0] = 0;
78867a12041SHong Zhang   nsend        = 0;
78915a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
79015a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
79115a3b8e2SHong Zhang     if (len_s[proc]) {
792445158ffSHong Zhang       nsend++;
79315a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
79415a3b8e2SHong Zhang       len         += len_si[proc];
79515a3b8e2SHong Zhang     }
79615a3b8e2SHong Zhang   }
79715a3b8e2SHong Zhang 
79815a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
799445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
800445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
80115a3b8e2SHong Zhang 
80215a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
80315a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
804445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
805445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
80615a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
80715a3b8e2SHong Zhang     if (!len_s[proc]) continue;
80815a3b8e2SHong Zhang     i    = owners_co[proc];
80915a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
81015a3b8e2SHong Zhang     k++;
81115a3b8e2SHong Zhang   }
81215a3b8e2SHong Zhang 
813681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
814681d504bSHong Zhang   /* ---------------------------------------- */
815681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
816681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
817681d504bSHong Zhang 
818681d504bSHong Zhang   /* receives coj are complete */
819445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
820445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
82115a3b8e2SHong Zhang   }
82215a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
823445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
82415a3b8e2SHong Zhang 
82578d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
82678d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
82778d30b94SHong Zhang     Jptr = buf_rj[k];
82878d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
82978d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
83078d30b94SHong Zhang     }
83178d30b94SHong Zhang   }
83278d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
83378d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
8349a6dcab7SHong Zhang 
83515a3b8e2SHong Zhang   /* (4) send and recv coi */
83615a3b8e2SHong Zhang   /*-----------------------*/
83715a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
838445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
83915a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
84015a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
84115a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
84215a3b8e2SHong Zhang     if (!len_s[proc]) continue;
84315a3b8e2SHong Zhang     /* form outgoing message for i-structure:
84415a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
84515a3b8e2SHong Zhang                [1:nrows]:           row index (global)
84615a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
84715a3b8e2SHong Zhang     */
84815a3b8e2SHong Zhang     /*-------------------------------------------*/
84915a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
85015a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
85115a3b8e2SHong Zhang     buf_si[0]   = nrows;
85215a3b8e2SHong Zhang     buf_si_i[0] = 0;
85315a3b8e2SHong Zhang     nrows       = 0;
85415a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
85515a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
85615a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
85715a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
85815a3b8e2SHong Zhang       nrows++;
85915a3b8e2SHong Zhang     }
86015a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
86115a3b8e2SHong Zhang     k++;
86215a3b8e2SHong Zhang     buf_si += len_si[proc];
86315a3b8e2SHong Zhang   }
864681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
865445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
86615a3b8e2SHong Zhang   }
86715a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
868445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
86915a3b8e2SHong Zhang 
8708cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
87115a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
87215a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
87315a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
874a4ffb656SHong Zhang 
87567a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
87667a12041SHong Zhang   /* ------------------------------------------ */
87778d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
87878d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
87915a3b8e2SHong Zhang   current_space = free_space;
88015a3b8e2SHong Zhang 
881445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
882445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
88315a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
88415a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
88515a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
88615a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
88715a3b8e2SHong Zhang   }
888445158ffSHong Zhang 
88978d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
89078d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
89115a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
89267a12041SHong Zhang     /* add C_loc into Cmpi */
89367a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
89467a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
89567a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
89615a3b8e2SHong Zhang 
89715a3b8e2SHong Zhang     /* add received col data into lnk */
898445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
89915a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
90015a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
90115a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
90215a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
90315a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
90415a3b8e2SHong Zhang       }
90515a3b8e2SHong Zhang     }
906d0e9a2f7SHong Zhang     nzi = lnk[0];
9078cb82516SHong Zhang 
90815a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
909d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
910d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
91115a3b8e2SHong Zhang   }
91215a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
91315a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
914445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
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;
936cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
9378cb82516SHong Zhang   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
9382259aa2eSHong Zhang 
9391a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9401a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9411a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
942aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
943cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
9441a47ec54SHong Zhang   *C                     = Cmpi;
9450d3441aeSHong Zhang   PetscFunctionReturn(0);
9460d3441aeSHong Zhang }
9470d3441aeSHong Zhang 
948*904d1e70Sandi selinger 
949aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
9500d3441aeSHong Zhang {
9510d3441aeSHong Zhang   PetscErrorCode    ierr;
9520d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
9530d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
954*904d1e70Sandi selinger   Mat_SeqAIJ        *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
9552dd9e643SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
9560d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
9579ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
9580d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
9598cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
9608cb82516SHong Zhang   PetscScalar       *apa;
9610d3441aeSHong Zhang   const PetscInt    *cols;
9620d3441aeSHong Zhang   const PetscScalar *vals;
9630d3441aeSHong Zhang 
9640d3441aeSHong Zhang   PetscFunctionBegin;
9650d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
966e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
967748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
9680d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
9690d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
970748c7196SHong Zhang   }
9710d3441aeSHong Zhang 
972e497d3c8SHong Zhang   /* 2) get AP_loc */
9730d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
9748cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
9750d3441aeSHong Zhang 
976e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
9770d3441aeSHong Zhang   /*-----------------------------------------------------*/
978748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
979748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
9800d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
9810d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
982e497d3c8SHong Zhang   }
9830d3441aeSHong Zhang 
9848cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
9858cb82516SHong Zhang   /* ---------------------------------------------- */
9860d3441aeSHong Zhang   /* get data from symbolic products */
9878cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
9882dd9e643SHong Zhang   if (ptap->P_oth) {
9898cb82516SHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
9902dd9e643SHong Zhang   }
9918cb82516SHong Zhang   apa   = ptap->apa;
992681d504bSHong Zhang   api   = ap->i;
993681d504bSHong Zhang   apj   = ap->j;
994e497d3c8SHong Zhang   for (i=0; i<am; i++) {
995445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
996e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
997e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
998e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
999e497d3c8SHong Zhang       col = apj[j+api[i]];
1000e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
1001e497d3c8SHong Zhang       apa[col] = 0.0;
10020d3441aeSHong Zhang     }
1003aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
10040d3441aeSHong Zhang   }
10050d3441aeSHong Zhang 
10068cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1007445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1008445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
10099ce11a7cSHong Zhang   C_loc = ptap->C_loc;
10109ce11a7cSHong Zhang   C_oth = ptap->C_oth;
10110d3441aeSHong Zhang 
10120d3441aeSHong Zhang   /* add C_loc and Co to to C */
10130d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
10140d3441aeSHong Zhang 
10150d3441aeSHong Zhang   /* C_loc -> C */
10160d3441aeSHong Zhang   cm    = C_loc->rmap->N;
10179ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
10188cb82516SHong Zhang   cols = c_seq->j;
10198cb82516SHong Zhang   vals = c_seq->a;
1020*904d1e70Sandi selinger 
1021*904d1e70Sandi selinger 
1022*904d1e70Sandi selinger   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assemled is PETSC_FALSE and */
1023*904d1e70Sandi selinger   /* when there are no off-processor parts.  */
1024*904d1e70Sandi selinger   if (C->assembled) {
1025*904d1e70Sandi selinger     C->was_assembled = PETSC_TRUE;
1026*904d1e70Sandi selinger     C->assembled     = PETSC_FALSE;
1027*904d1e70Sandi selinger   }
1028*904d1e70Sandi selinger   if (C->was_assembled) {
10290d3441aeSHong Zhang     for (i=0; i<cm; i++) {
10309ce11a7cSHong Zhang       ncols = c_seq->i[i+1] - c_seq->i[i];
10310d3441aeSHong Zhang       row = rstart + i;
1032*904d1e70Sandi selinger       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10338cb82516SHong Zhang       cols += ncols; vals += ncols;
10340d3441aeSHong Zhang     }
1035*904d1e70Sandi selinger   } else {
1036*904d1e70Sandi selinger     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a,cd->i,co->i);CHKERRQ(ierr);
1037*904d1e70Sandi selinger   }
10380d3441aeSHong Zhang 
10390d3441aeSHong Zhang   /* Co -> C, off-processor part */
10409ce11a7cSHong Zhang   cm = C_oth->rmap->N;
10419ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
10428cb82516SHong Zhang   cols = c_seq->j;
10438cb82516SHong Zhang   vals = c_seq->a;
10449ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
10459ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
10460d3441aeSHong Zhang     row = p->garray[i];
10470d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10488cb82516SHong Zhang     cols += ncols; vals += ncols;
10490d3441aeSHong Zhang   }
1050*904d1e70Sandi selinger 
10510d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10520d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10530d3441aeSHong Zhang 
1054748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
10550d3441aeSHong Zhang   PetscFunctionReturn(0);
10560d3441aeSHong Zhang }
1057