xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision b92f168fad6a7743cff583cc2bd4ec76d8b11400)
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) {
195419ecdd9Sandi selinger     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
196419ecdd9Sandi selinger     ierr = MatTranspose(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;
244edeac6deSandi selinger 
245e9ede7d0Sandi selinger   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
246edeac6deSandi selinger   /* when there are no off-processor parts.  */
2471de21080Sandi selinger   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
2481de21080Sandi selinger   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
2491de21080Sandi selinger   /* a table, and other, more complex stuff has to be done. */
250edeac6deSandi selinger   if (C->assembled) {
251edeac6deSandi selinger     C->was_assembled = PETSC_TRUE;
252edeac6deSandi selinger     C->assembled     = PETSC_FALSE;
253edeac6deSandi selinger   }
254edeac6deSandi selinger   if (C->was_assembled) {
255cf742fccSHong Zhang     for (i=0; i<cm; i++) {
256cf742fccSHong Zhang       ncols = c_seq->i[i+1] - c_seq->i[i];
257cf742fccSHong Zhang       row = rstart + i;
258edeac6deSandi selinger       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
259cf742fccSHong Zhang       cols += ncols; vals += ncols;
260cf742fccSHong Zhang     }
261edeac6deSandi selinger   } else {
262e9ede7d0Sandi selinger     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
263edeac6deSandi selinger   }
264cf742fccSHong Zhang 
265cf742fccSHong Zhang   /* Co -> C, off-processor part */
266cf742fccSHong Zhang   cm = C_oth->rmap->N;
267cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
268cf742fccSHong Zhang   cols = c_seq->j;
269cf742fccSHong Zhang   vals = c_seq->a;
270cf742fccSHong Zhang   for (i=0; i<cm; i++) {
271cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
272cf742fccSHong Zhang     row = p->garray[i];
273cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
274cf742fccSHong Zhang     cols += ncols; vals += ncols;
275cf742fccSHong Zhang   }
276cf742fccSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
277cf742fccSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
278cf742fccSHong Zhang 
279cf742fccSHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
280cf742fccSHong Zhang   PetscFunctionReturn(0);
281cf742fccSHong Zhang }
282cf742fccSHong Zhang 
283e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
2840d3441aeSHong Zhang {
2850d3441aeSHong Zhang   PetscErrorCode      ierr;
2860d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
287681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
2882259aa2eSHong Zhang   MPI_Comm            comm;
2892259aa2eSHong Zhang   PetscMPIInt         size,rank;
29076db11f6SHong Zhang   Mat                 Cmpi,P_loc,P_oth;
29115a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
292d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
293d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
294f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
29515a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
296d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
29715a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
29815a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
29915a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
300445158ffSHong Zhang   PetscLayout         rowmap;
301445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
302445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
303b4e8d1b6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi;
30452f7967eSHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
30567a12041SHong Zhang   PetscScalar         *apv;
30678d30b94SHong Zhang   PetscTable          ta;
307*b92f168fSBarry Smith   MatType             mtype;
308aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
3098cb82516SHong Zhang   PetscReal           apfill;
310aa690a28SHong Zhang #endif
31167a12041SHong Zhang 
31267a12041SHong Zhang   PetscFunctionBegin;
31367a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
31467a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
31567a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
316ae5f0867Sstefano_zampini 
31752f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
31852f7967eSHong Zhang 
319ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
320ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
321*b92f168fSBarry Smith   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
322*b92f168fSBarry Smith   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
323ae5f0867Sstefano_zampini 
324e953e47cSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
325e953e47cSHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
326e953e47cSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
327cf97cf3cSHong Zhang   ptap->algType = 0;
328e953e47cSHong Zhang 
329e953e47cSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
33076db11f6SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
331e953e47cSHong Zhang   /* get P_loc by taking all local rows of P */
33276db11f6SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
33376db11f6SHong Zhang 
33476db11f6SHong Zhang   ptap->P_loc = P_loc;
33576db11f6SHong Zhang   ptap->P_oth = P_oth;
336e953e47cSHong Zhang 
337e953e47cSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
338e953e47cSHong Zhang   /* --------------------------------- */
339419ecdd9Sandi selinger   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
340419ecdd9Sandi selinger   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
341e953e47cSHong Zhang 
342e953e47cSHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
343e953e47cSHong Zhang   /* ----------------------------------------------------------------- */
34476db11f6SHong Zhang   p_loc  = (Mat_SeqAIJ*)P_loc->data;
34552f7967eSHong Zhang   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
346e953e47cSHong Zhang 
347e953e47cSHong Zhang   /* create and initialize a linked list */
348e953e47cSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
34976db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
35076db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
351e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
352e953e47cSHong Zhang 
35376db11f6SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
354e953e47cSHong Zhang 
355e953e47cSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
35652f7967eSHong Zhang   if (ao) {
357e953e47cSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
35852f7967eSHong Zhang   } else {
35952f7967eSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
36052f7967eSHong Zhang   }
361e953e47cSHong Zhang   current_space = free_space;
362e953e47cSHong Zhang   nspacedouble  = 0;
363e953e47cSHong Zhang 
364e953e47cSHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
365e953e47cSHong Zhang   api[0] = 0;
366e953e47cSHong Zhang   for (i=0; i<am; i++) {
367e953e47cSHong Zhang     /* diagonal portion: Ad[i,:]*P */
368e953e47cSHong Zhang     ai = ad->i; pi = p_loc->i;
369e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
370e953e47cSHong Zhang     aj  = ad->j + ai[i];
371e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
372e953e47cSHong Zhang       row  = aj[j];
373e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
374e953e47cSHong Zhang       Jptr = p_loc->j + pi[row];
375e953e47cSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
37676db11f6SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
377e953e47cSHong Zhang     }
378e953e47cSHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
37952f7967eSHong Zhang     if (ao) {
380e953e47cSHong Zhang       ai = ao->i; pi = p_oth->i;
381e953e47cSHong Zhang       nzi = ai[i+1] - ai[i];
382e953e47cSHong Zhang       aj  = ao->j + ai[i];
383e953e47cSHong Zhang       for (j=0; j<nzi; j++) {
384e953e47cSHong Zhang         row  = aj[j];
385e953e47cSHong Zhang         pnz  = pi[row+1] - pi[row];
386e953e47cSHong Zhang         Jptr = p_oth->j + pi[row];
38776db11f6SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
388e953e47cSHong Zhang       }
38952f7967eSHong Zhang     }
390e953e47cSHong Zhang     apnz     = lnk[0];
391e953e47cSHong Zhang     api[i+1] = api[i] + apnz;
392e953e47cSHong Zhang 
393e953e47cSHong Zhang     /* if free space is not available, double the total space in the list */
394e953e47cSHong Zhang     if (current_space->local_remaining<apnz) {
395e953e47cSHong Zhang       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
396e953e47cSHong Zhang       nspacedouble++;
397e953e47cSHong Zhang     }
398e953e47cSHong Zhang 
399e953e47cSHong Zhang     /* Copy data into free space, then initialize lnk */
40076db11f6SHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
401e953e47cSHong Zhang 
402e953e47cSHong Zhang     current_space->array           += apnz;
403e953e47cSHong Zhang     current_space->local_used      += apnz;
404e953e47cSHong Zhang     current_space->local_remaining -= apnz;
405e953e47cSHong Zhang   }
406e953e47cSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
407e953e47cSHong Zhang   /* destroy list of free space and other temporary array(s) */
408e953e47cSHong Zhang   ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
409e953e47cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
41076db11f6SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
411e953e47cSHong Zhang 
412e953e47cSHong Zhang   /* Create AP_loc for reuse */
413e953e47cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
414e953e47cSHong Zhang 
415e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
41652f7967eSHong Zhang   if (ao) {
417e953e47cSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
41852f7967eSHong Zhang   } else {
41952f7967eSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
42052f7967eSHong Zhang   }
421e953e47cSHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
422e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
423e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
424e953e47cSHong Zhang 
425e953e47cSHong Zhang   if (api[am]) {
426b11c1ec8SHong 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);
427e953e47cSHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
428e953e47cSHong Zhang   } else {
429b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
430e953e47cSHong Zhang   }
431e953e47cSHong Zhang #endif
432e953e47cSHong Zhang 
433e953e47cSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
434e953e47cSHong Zhang   /* ------------------------------------ */
435e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
436e953e47cSHong Zhang 
437e953e47cSHong Zhang   /* (3) send coj of C_oth to other processors  */
438e953e47cSHong Zhang   /* ------------------------------------------ */
439e953e47cSHong Zhang   /* determine row ownership */
440e953e47cSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
441e953e47cSHong Zhang   rowmap->n  = pn;
442e953e47cSHong Zhang   rowmap->bs = 1;
443e953e47cSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
444e953e47cSHong Zhang   owners = rowmap->range;
445e953e47cSHong Zhang 
446e953e47cSHong Zhang   /* determine the number of messages to send, their lengths */
447e953e47cSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
448e953e47cSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
449e953e47cSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
450e953e47cSHong Zhang 
451e953e47cSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
452e953e47cSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
453e953e47cSHong Zhang   con   = ptap->C_oth->rmap->n;
454e953e47cSHong Zhang   proc  = 0;
455e953e47cSHong Zhang   for (i=0; i<con; i++) {
456e953e47cSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
457e953e47cSHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
458e953e47cSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
459e953e47cSHong Zhang   }
460e953e47cSHong Zhang 
461e953e47cSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
462e953e47cSHong Zhang   owners_co[0] = 0;
463e953e47cSHong Zhang   nsend        = 0;
464e953e47cSHong Zhang   for (proc=0; proc<size; proc++) {
465e953e47cSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
466e953e47cSHong Zhang     if (len_s[proc]) {
467e953e47cSHong Zhang       nsend++;
468e953e47cSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
469e953e47cSHong Zhang       len         += len_si[proc];
470e953e47cSHong Zhang     }
471e953e47cSHong Zhang   }
472e953e47cSHong Zhang 
473e953e47cSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
474e953e47cSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
475e953e47cSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
476e953e47cSHong Zhang 
477e953e47cSHong Zhang   /* post the Irecv and Isend of coj */
478e953e47cSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
479e953e47cSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
480e953e47cSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
481e953e47cSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
482e953e47cSHong Zhang     if (!len_s[proc]) continue;
483e953e47cSHong Zhang     i    = owners_co[proc];
484e953e47cSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
485e953e47cSHong Zhang     k++;
486e953e47cSHong Zhang   }
487e953e47cSHong Zhang 
488e953e47cSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
489e953e47cSHong Zhang   /* ---------------------------------------- */
490e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
491e953e47cSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
492e953e47cSHong Zhang 
493e953e47cSHong Zhang   /* receives coj are complete */
494e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
495e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
496e953e47cSHong Zhang   }
497e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
498e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
499e953e47cSHong Zhang 
500e953e47cSHong Zhang   /* add received column indices into ta to update Crmax */
501e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
502e953e47cSHong Zhang     Jptr = buf_rj[k];
503e953e47cSHong Zhang     for (j=0; j<len_r[k]; j++) {
504e953e47cSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
505e953e47cSHong Zhang     }
506e953e47cSHong Zhang   }
507e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
508e953e47cSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
509e953e47cSHong Zhang 
510e953e47cSHong Zhang   /* (4) send and recv coi */
511e953e47cSHong Zhang   /*-----------------------*/
512e953e47cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
513e953e47cSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
514e953e47cSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
515e953e47cSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
516e953e47cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
517e953e47cSHong Zhang     if (!len_s[proc]) continue;
518e953e47cSHong Zhang     /* form outgoing message for i-structure:
519e953e47cSHong Zhang          buf_si[0]:                 nrows to be sent
520e953e47cSHong Zhang                [1:nrows]:           row index (global)
521e953e47cSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
522e953e47cSHong Zhang     */
523e953e47cSHong Zhang     /*-------------------------------------------*/
524e953e47cSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
525e953e47cSHong Zhang     buf_si_i    = buf_si + nrows+1;
526e953e47cSHong Zhang     buf_si[0]   = nrows;
527e953e47cSHong Zhang     buf_si_i[0] = 0;
528e953e47cSHong Zhang     nrows       = 0;
529e953e47cSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
530e953e47cSHong Zhang       nzi = coi[i+1] - coi[i];
531e953e47cSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
532e953e47cSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
533e953e47cSHong Zhang       nrows++;
534e953e47cSHong Zhang     }
535e953e47cSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
536e953e47cSHong Zhang     k++;
537e953e47cSHong Zhang     buf_si += len_si[proc];
538e953e47cSHong Zhang   }
539e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
540e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
541e953e47cSHong Zhang   }
542e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
543e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
544e953e47cSHong Zhang 
545e953e47cSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
546e953e47cSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
547e953e47cSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
548e953e47cSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
549b4e8d1b6SHong Zhang 
550e953e47cSHong Zhang   /* (5) compute the local portion of Cmpi      */
551e953e47cSHong Zhang   /* ------------------------------------------ */
552e953e47cSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
553e953e47cSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
554e953e47cSHong Zhang   current_space = free_space;
555e953e47cSHong Zhang 
556e953e47cSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
557e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {
558e953e47cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
559e953e47cSHong Zhang     nrows       = *buf_ri_k[k];
560e953e47cSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
561e953e47cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
562e953e47cSHong Zhang   }
563e953e47cSHong Zhang 
564e953e47cSHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
565cf742fccSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
566e953e47cSHong Zhang   for (i=0; i<pn; i++) {
567e953e47cSHong Zhang     /* add C_loc into Cmpi */
568e953e47cSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
569e953e47cSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
570cf742fccSHong Zhang     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
571e953e47cSHong Zhang 
572e953e47cSHong Zhang     /* add received col data into lnk */
573e953e47cSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
574e953e47cSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
575e953e47cSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
576e953e47cSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
577cf742fccSHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
578e953e47cSHong Zhang         nextrow[k]++; nextci[k]++;
579e953e47cSHong Zhang       }
580e953e47cSHong Zhang     }
581e953e47cSHong Zhang     nzi = lnk[0];
582e953e47cSHong Zhang 
583e953e47cSHong Zhang     /* copy data into free space, then initialize lnk */
584cf742fccSHong Zhang     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
585e953e47cSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
586e953e47cSHong Zhang   }
587e953e47cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
588cf742fccSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
589e953e47cSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
590e953e47cSHong Zhang 
591e953e47cSHong Zhang   /* local sizes and preallocation */
592e953e47cSHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
593e953e47cSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
594e953e47cSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
595e953e47cSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
596e953e47cSHong Zhang 
597e953e47cSHong Zhang   /* members in merge */
598e953e47cSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
599e953e47cSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
600e953e47cSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
601e953e47cSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
602e953e47cSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
603e953e47cSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
604e953e47cSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
605e953e47cSHong Zhang 
606e953e47cSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
607e953e47cSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
608e953e47cSHong Zhang   c->ptap         = ptap;
609e953e47cSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
610e953e47cSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
611cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
612e953e47cSHong Zhang 
613e953e47cSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
614e953e47cSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
615a4ffb656SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
616e953e47cSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
617e953e47cSHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
618cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
619e953e47cSHong Zhang   *C                     = Cmpi;
620e953e47cSHong Zhang   PetscFunctionReturn(0);
621e953e47cSHong Zhang }
622e953e47cSHong Zhang 
623e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
624e953e47cSHong Zhang {
625e953e47cSHong Zhang   PetscErrorCode      ierr;
626e953e47cSHong Zhang   Mat_PtAPMPI         *ptap;
627e953e47cSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
628e953e47cSHong Zhang   MPI_Comm            comm;
629e953e47cSHong Zhang   PetscMPIInt         size,rank;
630e953e47cSHong Zhang   Mat                 Cmpi;
631e953e47cSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
632e953e47cSHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
633e953e47cSHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
634e953e47cSHong Zhang   PetscBT             lnkbt;
635e953e47cSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
636e953e47cSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
637e953e47cSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
638e953e47cSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
639e953e47cSHong Zhang   MPI_Request         *swaits,*rwaits;
640e953e47cSHong Zhang   MPI_Status          *sstatus,rstatus;
641e953e47cSHong Zhang   PetscLayout         rowmap;
642e953e47cSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
643e953e47cSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
644e953e47cSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
645ec07b8f8SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
646e953e47cSHong Zhang   PetscScalar         *apv;
647e953e47cSHong Zhang   PetscTable          ta;
648*b92f168fSBarry Smith   MatType             mtype;
649e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
650e953e47cSHong Zhang   PetscReal           apfill;
651e953e47cSHong Zhang #endif
652e953e47cSHong Zhang 
653e953e47cSHong Zhang   PetscFunctionBegin;
654b4e8d1b6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
655e953e47cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
656e953e47cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
657e953e47cSHong Zhang 
65852f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
659ec07b8f8SHong Zhang 
660e953e47cSHong Zhang   /* create symbolic parallel matrix Cmpi */
661e953e47cSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
662*b92f168fSBarry Smith   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
663*b92f168fSBarry Smith   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
664e953e47cSHong Zhang 
665e953e47cSHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
666e953e47cSHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
667e953e47cSHong Zhang 
66815a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
66915a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
67015a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
671a4ffb656SHong Zhang   ptap->algType = 1;
67215a3b8e2SHong Zhang 
67315a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
67415a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
67515a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
67615a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
67715a3b8e2SHong Zhang 
67867a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
67967a12041SHong Zhang   /* --------------------------------- */
680419ecdd9Sandi selinger   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
681419ecdd9Sandi selinger   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
68215a3b8e2SHong Zhang 
68367a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
68467a12041SHong Zhang   /* ----------------------------------------------------------------- */
68567a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
68652f7967eSHong Zhang   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
687445158ffSHong Zhang 
6889a6dcab7SHong Zhang   /* create and initialize a linked list */
68945d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
6904b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
6914b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
69278d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
693d0e9a2f7SHong 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); */
69478d30b94SHong Zhang 
69578d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
696445158ffSHong Zhang 
6978cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
698ec07b8f8SHong Zhang   if (ao) {
699f91af8c7SBarry Smith     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
700ec07b8f8SHong Zhang   } else {
701ec07b8f8SHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
702ec07b8f8SHong Zhang   }
703445158ffSHong Zhang   current_space = free_space;
70467a12041SHong Zhang   nspacedouble  = 0;
70567a12041SHong Zhang 
70667a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
70767a12041SHong Zhang   api[0] = 0;
708445158ffSHong Zhang   for (i=0; i<am; i++) {
70967a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
71067a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
71167a12041SHong Zhang     nzi = ai[i+1] - ai[i];
71267a12041SHong Zhang     aj  = ad->j + ai[i];
713445158ffSHong Zhang     for (j=0; j<nzi; j++) {
714445158ffSHong Zhang       row  = aj[j];
71567a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
71667a12041SHong Zhang       Jptr = p_loc->j + pi[row];
717445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
718445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
719445158ffSHong Zhang     }
72067a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
721ec07b8f8SHong Zhang     if (ao) {
72267a12041SHong Zhang       ai = ao->i; pi = p_oth->i;
72367a12041SHong Zhang       nzi = ai[i+1] - ai[i];
72467a12041SHong Zhang       aj  = ao->j + ai[i];
725445158ffSHong Zhang       for (j=0; j<nzi; j++) {
726445158ffSHong Zhang         row  = aj[j];
72767a12041SHong Zhang         pnz  = pi[row+1] - pi[row];
72867a12041SHong Zhang         Jptr = p_oth->j + pi[row];
729445158ffSHong Zhang         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
730445158ffSHong Zhang       }
731ec07b8f8SHong Zhang     }
732445158ffSHong Zhang     apnz     = lnk[0];
733445158ffSHong Zhang     api[i+1] = api[i] + apnz;
734445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
735445158ffSHong Zhang 
736445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
737445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
738f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
739445158ffSHong Zhang       nspacedouble++;
740445158ffSHong Zhang     }
741445158ffSHong Zhang 
742445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
743445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
744445158ffSHong Zhang 
745445158ffSHong Zhang     current_space->array           += apnz;
746445158ffSHong Zhang     current_space->local_used      += apnz;
747445158ffSHong Zhang     current_space->local_remaining -= apnz;
748445158ffSHong Zhang   }
749681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
750445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
751445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
752445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
7539a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
7549a6dcab7SHong Zhang 
755aa690a28SHong Zhang   /* Create AP_loc for reuse */
756445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
757aa690a28SHong Zhang 
758aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
759ec07b8f8SHong Zhang   if (ao) {
760aa690a28SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
761ec07b8f8SHong Zhang   } else {
762ec07b8f8SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
763ec07b8f8SHong Zhang   }
764aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
765aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
766aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
767aa690a28SHong Zhang 
768aa690a28SHong Zhang   if (api[am]) {
769b11c1ec8SHong 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);
770aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
771aa690a28SHong Zhang   } else {
772b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
773aa690a28SHong Zhang   }
774aa690a28SHong Zhang #endif
775aa690a28SHong Zhang 
776681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
777681d504bSHong Zhang   /* ------------------------------------ */
77867a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
77915a3b8e2SHong Zhang 
78067a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
78167a12041SHong Zhang   /* ------------------------------------------ */
78215a3b8e2SHong Zhang   /* determine row ownership */
783445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
784445158ffSHong Zhang   rowmap->n  = pn;
785445158ffSHong Zhang   rowmap->bs = 1;
786445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
787445158ffSHong Zhang   owners = rowmap->range;
78815a3b8e2SHong Zhang 
78915a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
7908cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
7918cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
79215a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
79315a3b8e2SHong Zhang 
79467a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
79567a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
79667a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
79715a3b8e2SHong Zhang   proc  = 0;
79867a12041SHong Zhang   for (i=0; i<con; i++) {
79915a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
80015a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
80115a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
80215a3b8e2SHong Zhang   }
80315a3b8e2SHong Zhang 
80415a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
80515a3b8e2SHong Zhang   owners_co[0] = 0;
80667a12041SHong Zhang   nsend        = 0;
80715a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
80815a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
80915a3b8e2SHong Zhang     if (len_s[proc]) {
810445158ffSHong Zhang       nsend++;
81115a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
81215a3b8e2SHong Zhang       len         += len_si[proc];
81315a3b8e2SHong Zhang     }
81415a3b8e2SHong Zhang   }
81515a3b8e2SHong Zhang 
81615a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
817445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
818445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
81915a3b8e2SHong Zhang 
82015a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
82115a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
822445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
823445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
82415a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
82515a3b8e2SHong Zhang     if (!len_s[proc]) continue;
82615a3b8e2SHong Zhang     i    = owners_co[proc];
82715a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
82815a3b8e2SHong Zhang     k++;
82915a3b8e2SHong Zhang   }
83015a3b8e2SHong Zhang 
831681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
832681d504bSHong Zhang   /* ---------------------------------------- */
833681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
834681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
835681d504bSHong Zhang 
836681d504bSHong Zhang   /* receives coj are complete */
837445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
838445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
83915a3b8e2SHong Zhang   }
84015a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
841445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
84215a3b8e2SHong Zhang 
84378d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
84478d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
84578d30b94SHong Zhang     Jptr = buf_rj[k];
84678d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
84778d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
84878d30b94SHong Zhang     }
84978d30b94SHong Zhang   }
85078d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
85178d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
8529a6dcab7SHong Zhang 
85315a3b8e2SHong Zhang   /* (4) send and recv coi */
85415a3b8e2SHong Zhang   /*-----------------------*/
85515a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
856445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
85715a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
85815a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
85915a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
86015a3b8e2SHong Zhang     if (!len_s[proc]) continue;
86115a3b8e2SHong Zhang     /* form outgoing message for i-structure:
86215a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
86315a3b8e2SHong Zhang                [1:nrows]:           row index (global)
86415a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
86515a3b8e2SHong Zhang     */
86615a3b8e2SHong Zhang     /*-------------------------------------------*/
86715a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
86815a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
86915a3b8e2SHong Zhang     buf_si[0]   = nrows;
87015a3b8e2SHong Zhang     buf_si_i[0] = 0;
87115a3b8e2SHong Zhang     nrows       = 0;
87215a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
87315a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
87415a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
87515a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
87615a3b8e2SHong Zhang       nrows++;
87715a3b8e2SHong Zhang     }
87815a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
87915a3b8e2SHong Zhang     k++;
88015a3b8e2SHong Zhang     buf_si += len_si[proc];
88115a3b8e2SHong Zhang   }
882681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
883445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
88415a3b8e2SHong Zhang   }
88515a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
886445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
88715a3b8e2SHong Zhang 
8888cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
88915a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
89015a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
89115a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
892a4ffb656SHong Zhang 
89367a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
89467a12041SHong Zhang   /* ------------------------------------------ */
89578d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
89678d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
89715a3b8e2SHong Zhang   current_space = free_space;
89815a3b8e2SHong Zhang 
899445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
900445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
90115a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
90215a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
90315a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
90415a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
90515a3b8e2SHong Zhang   }
906445158ffSHong Zhang 
90778d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
90878d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
90915a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
91067a12041SHong Zhang     /* add C_loc into Cmpi */
91167a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
91267a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
91367a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
91415a3b8e2SHong Zhang 
91515a3b8e2SHong Zhang     /* add received col data into lnk */
916445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
91715a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
91815a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
91915a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
92015a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
92115a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
92215a3b8e2SHong Zhang       }
92315a3b8e2SHong Zhang     }
924d0e9a2f7SHong Zhang     nzi = lnk[0];
9258cb82516SHong Zhang 
92615a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
927d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
928d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
92915a3b8e2SHong Zhang   }
93015a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
93115a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
932445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
93380bb4639SHong Zhang 
934ae5f0867Sstefano_zampini   /* local sizes and preallocation */
93515a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
93615a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
93715a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
93815a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
93915a3b8e2SHong Zhang 
940445158ffSHong Zhang   /* members in merge */
941445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
942445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
943445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
944445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
945445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
946445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
947445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
94815a3b8e2SHong Zhang 
94915a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
95015a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
95115a3b8e2SHong Zhang   c->ptap         = ptap;
9521a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
9531a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
954cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
9558cb82516SHong Zhang   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
9562259aa2eSHong Zhang 
9571a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9581a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9591a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
960aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
961cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
9621a47ec54SHong Zhang   *C                     = Cmpi;
9630d3441aeSHong Zhang   PetscFunctionReturn(0);
9640d3441aeSHong Zhang }
9650d3441aeSHong Zhang 
966904d1e70Sandi selinger 
967aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
9680d3441aeSHong Zhang {
9690d3441aeSHong Zhang   PetscErrorCode    ierr;
9700d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
9710d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
9722dd9e643SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
9730d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
9749ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
9750d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
9768cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
9778cb82516SHong Zhang   PetscScalar       *apa;
9780d3441aeSHong Zhang   const PetscInt    *cols;
9790d3441aeSHong Zhang   const PetscScalar *vals;
9800d3441aeSHong Zhang 
9810d3441aeSHong Zhang   PetscFunctionBegin;
9820d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
983e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
984748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
985419ecdd9Sandi selinger     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
986419ecdd9Sandi selinger     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
987748c7196SHong Zhang   }
9880d3441aeSHong Zhang 
989e497d3c8SHong Zhang   /* 2) get AP_loc */
9900d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
9918cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
9920d3441aeSHong Zhang 
993e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
9940d3441aeSHong Zhang   /*-----------------------------------------------------*/
995748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
996748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
9970d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
9980d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
999e497d3c8SHong Zhang   }
10000d3441aeSHong Zhang 
10018cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
10028cb82516SHong Zhang   /* ---------------------------------------------- */
10030d3441aeSHong Zhang   /* get data from symbolic products */
10048cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
10052dd9e643SHong Zhang   if (ptap->P_oth) {
10068cb82516SHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
10072dd9e643SHong Zhang   }
10088cb82516SHong Zhang   apa   = ptap->apa;
1009681d504bSHong Zhang   api   = ap->i;
1010681d504bSHong Zhang   apj   = ap->j;
1011e497d3c8SHong Zhang   for (i=0; i<am; i++) {
1012445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
1013e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
1014e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
1015e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
1016e497d3c8SHong Zhang       col = apj[j+api[i]];
1017e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
1018e497d3c8SHong Zhang       apa[col] = 0.0;
10190d3441aeSHong Zhang     }
1020aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
10210d3441aeSHong Zhang   }
10220d3441aeSHong Zhang 
10238cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1024445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1025445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
10269ce11a7cSHong Zhang   C_loc = ptap->C_loc;
10279ce11a7cSHong Zhang   C_oth = ptap->C_oth;
10280d3441aeSHong Zhang 
10290d3441aeSHong Zhang   /* add C_loc and Co to to C */
10300d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
10310d3441aeSHong Zhang 
10320d3441aeSHong Zhang   /* C_loc -> C */
10330d3441aeSHong Zhang   cm    = C_loc->rmap->N;
10349ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
10358cb82516SHong Zhang   cols = c_seq->j;
10368cb82516SHong Zhang   vals = c_seq->a;
1037904d1e70Sandi selinger 
1038904d1e70Sandi selinger 
1039e9ede7d0Sandi selinger   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
1040904d1e70Sandi selinger   /* when there are no off-processor parts.  */
10411de21080Sandi selinger   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
10421de21080Sandi selinger   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
10431de21080Sandi selinger   /* a table, and other, more complex stuff has to be done. */
1044904d1e70Sandi selinger   if (C->assembled) {
1045904d1e70Sandi selinger     C->was_assembled = PETSC_TRUE;
1046904d1e70Sandi selinger     C->assembled     = PETSC_FALSE;
1047904d1e70Sandi selinger   }
1048904d1e70Sandi selinger   if (C->was_assembled) {
10490d3441aeSHong Zhang     for (i=0; i<cm; i++) {
10509ce11a7cSHong Zhang       ncols = c_seq->i[i+1] - c_seq->i[i];
10510d3441aeSHong Zhang       row = rstart + i;
1052904d1e70Sandi selinger       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10538cb82516SHong Zhang       cols += ncols; vals += ncols;
10540d3441aeSHong Zhang     }
1055904d1e70Sandi selinger   } else {
1056e9ede7d0Sandi selinger     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
1057904d1e70Sandi selinger   }
10580d3441aeSHong Zhang 
10590d3441aeSHong Zhang   /* Co -> C, off-processor part */
10609ce11a7cSHong Zhang   cm = C_oth->rmap->N;
10619ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
10628cb82516SHong Zhang   cols = c_seq->j;
10638cb82516SHong Zhang   vals = c_seq->a;
10649ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
10659ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
10660d3441aeSHong Zhang     row = p->garray[i];
10670d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10688cb82516SHong Zhang     cols += ncols; vals += ncols;
10690d3441aeSHong Zhang   }
1070904d1e70Sandi selinger 
10710d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10720d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10730d3441aeSHong Zhang 
1074748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
10750d3441aeSHong Zhang   PetscFunctionReturn(0);
10760d3441aeSHong Zhang }
1077