xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision a9899c979f2dda621d6accdc0838ebb8823aa852)
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 
393697aca6SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAPMPI_private(Mat_PtAPMPI *ptap)
40cc6cb787SHong Zhang {
41cc6cb787SHong Zhang   PetscErrorCode      ierr;
423697aca6SHong Zhang   Mat_Merge_SeqsToMPI *merge=ptap->merge;
43cc6cb787SHong Zhang 
44cc6cb787SHong Zhang   PetscFunctionBegin;
45b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
46f8487c73SHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
47a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
48a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
49c5af039cSHong Zhang   ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
5005d62848SHong Zhang   ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
5105d62848SHong Zhang   ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
528403a639SHong Zhang   if (ptap->AP_loc) { /* used by alg_rap */
53681d504bSHong Zhang     Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
54681d504bSHong Zhang     ierr = PetscFree(ap->i);CHKERRQ(ierr);
55681d504bSHong Zhang     ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
560d3441aeSHong Zhang     ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
578403a639SHong Zhang   } else { /* used by alg_ptap */
588403a639SHong Zhang     ierr = PetscFree(ptap->api);CHKERRQ(ierr);
598403a639SHong Zhang     ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
60681d504bSHong Zhang   }
612259aa2eSHong Zhang   ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
622259aa2eSHong Zhang   ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
63414894bdSHong Zhang   if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
64a560ef98SHong Zhang 
658403a639SHong Zhang   if (merge) { /* used by alg_ptap */
66cc6cb787SHong Zhang     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
67cc6cb787SHong Zhang     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
68cc6cb787SHong Zhang     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
69cc6cb787SHong Zhang     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
70cc6cb787SHong Zhang     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
71c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
72cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
73c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
74cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
75445158ffSHong Zhang     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
76445158ffSHong Zhang     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
7705b42c5fSBarry Smith     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
786bf464f9SBarry Smith     ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
79f8487c73SHong Zhang     ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
80bf0cc555SLisandro Dalcin   }
813697aca6SHong Zhang   PetscFunctionReturn(0);
823697aca6SHong Zhang }
83a560ef98SHong Zhang 
843697aca6SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
853697aca6SHong Zhang {
863697aca6SHong Zhang   PetscErrorCode ierr;
873697aca6SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
883697aca6SHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
893697aca6SHong Zhang 
903697aca6SHong Zhang   PetscFunctionBegin;
913697aca6SHong Zhang   if (ptap) {
923697aca6SHong Zhang     if (A->reuse) {
933697aca6SHong Zhang       ierr = MatDestroy_MPIAIJ_PtAPMPI_private(a->ptap);CHKERRQ(ierr);
943697aca6SHong Zhang     }
95a560ef98SHong Zhang     ierr = ptap->destroy(A);CHKERRQ(ierr);
96f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
97c8b0795cSMark F. Adams   }
98cc6cb787SHong Zhang   PetscFunctionReturn(0);
99cc6cb787SHong Zhang }
100cc6cb787SHong Zhang 
101150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
10242c57489SHong Zhang {
10342c57489SHong Zhang   PetscErrorCode ierr;
104a4ffb656SHong Zhang   PetscBool      flg;
10567a12041SHong Zhang   MPI_Comm       comm;
106a4ffb656SHong Zhang #if !defined(PETSC_HAVE_HYPRE)
107a4ffb656SHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
108a4ffb656SHong Zhang   PetscInt            nalg=2;
109a4ffb656SHong Zhang #else
110a4ffb656SHong Zhang   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
111a4ffb656SHong Zhang   PetscInt            nalg=3;
112a4ffb656SHong Zhang #endif
113a4ffb656SHong Zhang   PetscInt            pN=P->cmap->N,alg=1; /* set default algorithm */
11442c57489SHong Zhang 
11542c57489SHong Zhang   PetscFunctionBegin;
11667a12041SHong Zhang   /* check if matrix local sizes are compatible */
11767a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1186c4ed002SBarry 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);
1196c4ed002SBarry 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);
12067a12041SHong Zhang 
121cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
122a4ffb656SHong Zhang     /* pick an algorithm */
123a4ffb656SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
124a4ffb656SHong Zhang     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
125a4ffb656SHong Zhang     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
126a4ffb656SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
127a4ffb656SHong Zhang 
128a4ffb656SHong Zhang     if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
129a4ffb656SHong Zhang       MatInfo     Ainfo,Pinfo;
130a4ffb656SHong Zhang       PetscInt    nz_local;
131a4ffb656SHong Zhang       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
132a4ffb656SHong Zhang 
133a4ffb656SHong Zhang       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
134a4ffb656SHong Zhang       ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
135a4ffb656SHong Zhang       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
136a4ffb656SHong Zhang 
137a4ffb656SHong Zhang       if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
138a4ffb656SHong Zhang       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
139a4ffb656SHong Zhang 
140a4ffb656SHong Zhang       if (alg_scalable) {
141a4ffb656SHong Zhang         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
1420d3441aeSHong Zhang       }
143a4ffb656SHong Zhang     }
144a4ffb656SHong Zhang 
145a4ffb656SHong Zhang     switch (alg) {
146a4ffb656SHong Zhang     case 1:
147a4ffb656SHong Zhang       /* do R=P^T locally, then C=R*A*P */
148a4ffb656SHong Zhang       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
149a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
1503ff4c91cSHong Zhang       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
151a4ffb656SHong Zhang       break;
152a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE)
153a4ffb656SHong Zhang     case 2:
154a4ffb656SHong Zhang       /* Use boomerAMGBuildCoarseOperator */
155a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
156a4ffb656SHong Zhang       PetscFunctionReturn(0);
157a4ffb656SHong Zhang       break;
158a4ffb656SHong Zhang #endif
159a4ffb656SHong Zhang     default:
160a4ffb656SHong Zhang       /* do R=P^T locally, then C=R*A*P */
161a4ffb656SHong Zhang       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
162a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
163a4ffb656SHong Zhang       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
164a4ffb656SHong Zhang       break;
165a4ffb656SHong Zhang     }
1667a7894deSKris Buschelman   }
1673ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1688403a639SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
1693ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
17042c57489SHong Zhang   PetscFunctionReturn(0);
17142c57489SHong Zhang }
17242c57489SHong Zhang 
173cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
174cf742fccSHong Zhang {
175cf742fccSHong Zhang   PetscErrorCode    ierr;
176cf742fccSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
177cf742fccSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
178cf742fccSHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
179cf742fccSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
180cf742fccSHong Zhang   Mat               AP_loc,C_loc,C_oth;
1815ca39647SHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz;
182cf742fccSHong Zhang   PetscScalar       *apa;
183cf742fccSHong Zhang   const PetscInt    *cols;
184cf742fccSHong Zhang   const PetscScalar *vals;
185cf742fccSHong Zhang 
186cf742fccSHong Zhang   PetscFunctionBegin;
187cf742fccSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
188cf742fccSHong Zhang 
189cf742fccSHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
190cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
191419ecdd9Sandi selinger     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
192419ecdd9Sandi selinger     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
193cf742fccSHong Zhang   }
194cf742fccSHong Zhang 
195cf742fccSHong Zhang   /* 2) get AP_loc */
196cf742fccSHong Zhang   AP_loc = ptap->AP_loc;
197cf742fccSHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
198cf742fccSHong Zhang 
199cf742fccSHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
200cf742fccSHong Zhang   /*-----------------------------------------------------*/
201cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
202cf742fccSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
203cf742fccSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
204cf742fccSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
205cf742fccSHong Zhang   }
206cf742fccSHong Zhang 
207cf742fccSHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
208cf742fccSHong Zhang   /* ---------------------------------------------- */
209cf742fccSHong Zhang   /* get data from symbolic products */
210cf742fccSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
21152f7967eSHong Zhang   if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
212c072d3e2SSatish Balay   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!");
21352f7967eSHong Zhang 
214cf742fccSHong Zhang   api   = ap->i;
215cf742fccSHong Zhang   apj   = ap->j;
216cf742fccSHong Zhang   for (i=0; i<am; i++) {
217cf742fccSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
218cf742fccSHong Zhang     apnz = api[i+1] - api[i];
219b4e8d1b6SHong Zhang     apa = ap->a + api[i];
220b4e8d1b6SHong Zhang     ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
221b4e8d1b6SHong Zhang     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
222cf742fccSHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
223cf742fccSHong Zhang   }
224cf742fccSHong Zhang 
225cf742fccSHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
226cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
227cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
228cf742fccSHong Zhang 
229cf742fccSHong Zhang   C_loc = ptap->C_loc;
230cf742fccSHong Zhang   C_oth = ptap->C_oth;
231cf742fccSHong Zhang 
232cf742fccSHong Zhang   /* add C_loc and Co to to C */
233cf742fccSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
234cf742fccSHong Zhang 
235cf742fccSHong Zhang   /* C_loc -> C */
236cf742fccSHong Zhang   cm    = C_loc->rmap->N;
237cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
238cf742fccSHong Zhang   cols = c_seq->j;
239cf742fccSHong Zhang   vals = c_seq->a;
240edeac6deSandi selinger 
241e9ede7d0Sandi selinger   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
242edeac6deSandi selinger   /* when there are no off-processor parts.  */
2431de21080Sandi selinger   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
2441de21080Sandi selinger   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
2451de21080Sandi selinger   /* a table, and other, more complex stuff has to be done. */
246edeac6deSandi selinger   if (C->assembled) {
247edeac6deSandi selinger     C->was_assembled = PETSC_TRUE;
248edeac6deSandi selinger     C->assembled     = PETSC_FALSE;
249edeac6deSandi selinger   }
250edeac6deSandi selinger   if (C->was_assembled) {
251cf742fccSHong Zhang     for (i=0; i<cm; i++) {
252cf742fccSHong Zhang       ncols = c_seq->i[i+1] - c_seq->i[i];
253cf742fccSHong Zhang       row = rstart + i;
254edeac6deSandi selinger       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
255cf742fccSHong Zhang       cols += ncols; vals += ncols;
256cf742fccSHong Zhang     }
257edeac6deSandi selinger   } else {
258e9ede7d0Sandi selinger     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
259edeac6deSandi selinger   }
260cf742fccSHong Zhang 
261cf742fccSHong Zhang   /* Co -> C, off-processor part */
262cf742fccSHong Zhang   cm = C_oth->rmap->N;
263cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
264cf742fccSHong Zhang   cols = c_seq->j;
265cf742fccSHong Zhang   vals = c_seq->a;
266cf742fccSHong Zhang   for (i=0; i<cm; i++) {
267cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
268cf742fccSHong Zhang     row = p->garray[i];
269cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
270cf742fccSHong Zhang     cols += ncols; vals += ncols;
271cf742fccSHong Zhang   }
272cf742fccSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
273cf742fccSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
274cf742fccSHong Zhang 
275cf742fccSHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
276cf742fccSHong Zhang   PetscFunctionReturn(0);
277cf742fccSHong Zhang }
278cf742fccSHong Zhang 
279e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
2800d3441aeSHong Zhang {
2810d3441aeSHong Zhang   PetscErrorCode      ierr;
2820d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
283681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
2842259aa2eSHong Zhang   MPI_Comm            comm;
2852259aa2eSHong Zhang   PetscMPIInt         size,rank;
28676db11f6SHong Zhang   Mat                 Cmpi,P_loc,P_oth;
28715a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
288d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
289d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
290f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
29115a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
292d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
29315a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
29415a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
29515a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
296445158ffSHong Zhang   PetscLayout         rowmap;
297445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
298445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
299b4e8d1b6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi;
30052f7967eSHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
30167a12041SHong Zhang   PetscScalar         *apv;
30278d30b94SHong Zhang   PetscTable          ta;
303b92f168fSBarry Smith   MatType             mtype;
304aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
3058cb82516SHong Zhang   PetscReal           apfill;
306aa690a28SHong Zhang #endif
30767a12041SHong Zhang 
30867a12041SHong Zhang   PetscFunctionBegin;
30967a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
31067a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
31167a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
312ae5f0867Sstefano_zampini 
31352f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
31452f7967eSHong Zhang 
315ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
316ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
317b92f168fSBarry Smith   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
318b92f168fSBarry Smith   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
319ae5f0867Sstefano_zampini 
320e953e47cSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
321e953e47cSHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
322e953e47cSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
323cf97cf3cSHong Zhang   ptap->algType = 0;
324e953e47cSHong Zhang 
325e953e47cSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
32676db11f6SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
327e953e47cSHong Zhang   /* get P_loc by taking all local rows of P */
32876db11f6SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
32976db11f6SHong Zhang 
33076db11f6SHong Zhang   ptap->P_loc = P_loc;
33176db11f6SHong Zhang   ptap->P_oth = P_oth;
332e953e47cSHong Zhang 
333e953e47cSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
334e953e47cSHong Zhang   /* --------------------------------- */
335419ecdd9Sandi selinger   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
336419ecdd9Sandi selinger   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
337e953e47cSHong Zhang 
338e953e47cSHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
339e953e47cSHong Zhang   /* ----------------------------------------------------------------- */
34076db11f6SHong Zhang   p_loc  = (Mat_SeqAIJ*)P_loc->data;
34152f7967eSHong Zhang   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
342e953e47cSHong Zhang 
343e953e47cSHong Zhang   /* create and initialize a linked list */
344e953e47cSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
34576db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
34676db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
347e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
348e953e47cSHong Zhang 
34976db11f6SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
350e953e47cSHong Zhang 
351e953e47cSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
35252f7967eSHong Zhang   if (ao) {
353e953e47cSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
35452f7967eSHong Zhang   } else {
35552f7967eSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
35652f7967eSHong Zhang   }
357e953e47cSHong Zhang   current_space = free_space;
358e953e47cSHong Zhang   nspacedouble  = 0;
359e953e47cSHong Zhang 
360e953e47cSHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
361e953e47cSHong Zhang   api[0] = 0;
362e953e47cSHong Zhang   for (i=0; i<am; i++) {
363e953e47cSHong Zhang     /* diagonal portion: Ad[i,:]*P */
364e953e47cSHong Zhang     ai = ad->i; pi = p_loc->i;
365e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
366e953e47cSHong Zhang     aj  = ad->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_loc->j + pi[row];
371e953e47cSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
37276db11f6SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
373e953e47cSHong Zhang     }
374e953e47cSHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
37552f7967eSHong Zhang     if (ao) {
376e953e47cSHong Zhang       ai = ao->i; pi = p_oth->i;
377e953e47cSHong Zhang       nzi = ai[i+1] - ai[i];
378e953e47cSHong Zhang       aj  = ao->j + ai[i];
379e953e47cSHong Zhang       for (j=0; j<nzi; j++) {
380e953e47cSHong Zhang         row  = aj[j];
381e953e47cSHong Zhang         pnz  = pi[row+1] - pi[row];
382e953e47cSHong Zhang         Jptr = p_oth->j + pi[row];
38376db11f6SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
384e953e47cSHong Zhang       }
38552f7967eSHong Zhang     }
386e953e47cSHong Zhang     apnz     = lnk[0];
387e953e47cSHong Zhang     api[i+1] = api[i] + apnz;
388e953e47cSHong Zhang 
389e953e47cSHong Zhang     /* if free space is not available, double the total space in the list */
390e953e47cSHong Zhang     if (current_space->local_remaining<apnz) {
391e953e47cSHong Zhang       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
392e953e47cSHong Zhang       nspacedouble++;
393e953e47cSHong Zhang     }
394e953e47cSHong Zhang 
395e953e47cSHong Zhang     /* Copy data into free space, then initialize lnk */
39676db11f6SHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
397e953e47cSHong Zhang 
398e953e47cSHong Zhang     current_space->array           += apnz;
399e953e47cSHong Zhang     current_space->local_used      += apnz;
400e953e47cSHong Zhang     current_space->local_remaining -= apnz;
401e953e47cSHong Zhang   }
402e953e47cSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
403e953e47cSHong Zhang   /* destroy list of free space and other temporary array(s) */
404e953e47cSHong Zhang   ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
405e953e47cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
40676db11f6SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
407e953e47cSHong Zhang 
408e953e47cSHong Zhang   /* Create AP_loc for reuse */
409e953e47cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
410e953e47cSHong Zhang 
411e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
41252f7967eSHong Zhang   if (ao) {
413e953e47cSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
41452f7967eSHong Zhang   } else {
41552f7967eSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
41652f7967eSHong Zhang   }
417e953e47cSHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
418e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
419e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
420e953e47cSHong Zhang 
421e953e47cSHong Zhang   if (api[am]) {
422b11c1ec8SHong 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);
423e953e47cSHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
424e953e47cSHong Zhang   } else {
425b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
426e953e47cSHong Zhang   }
427e953e47cSHong Zhang #endif
428e953e47cSHong Zhang 
429e953e47cSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
430e953e47cSHong Zhang   /* ------------------------------------ */
431e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
432e953e47cSHong Zhang 
433e953e47cSHong Zhang   /* (3) send coj of C_oth to other processors  */
434e953e47cSHong Zhang   /* ------------------------------------------ */
435e953e47cSHong Zhang   /* determine row ownership */
436e953e47cSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
437e953e47cSHong Zhang   rowmap->n  = pn;
438e953e47cSHong Zhang   rowmap->bs = 1;
439e953e47cSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
440e953e47cSHong Zhang   owners = rowmap->range;
441e953e47cSHong Zhang 
442e953e47cSHong Zhang   /* determine the number of messages to send, their lengths */
443e953e47cSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
444e953e47cSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
445e953e47cSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
446e953e47cSHong Zhang 
447e953e47cSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
448e953e47cSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
449e953e47cSHong Zhang   con   = ptap->C_oth->rmap->n;
450e953e47cSHong Zhang   proc  = 0;
451e953e47cSHong Zhang   for (i=0; i<con; i++) {
452e953e47cSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
453e953e47cSHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
454e953e47cSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
455e953e47cSHong Zhang   }
456e953e47cSHong Zhang 
457e953e47cSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
458e953e47cSHong Zhang   owners_co[0] = 0;
459e953e47cSHong Zhang   nsend        = 0;
460e953e47cSHong Zhang   for (proc=0; proc<size; proc++) {
461e953e47cSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
462e953e47cSHong Zhang     if (len_s[proc]) {
463e953e47cSHong Zhang       nsend++;
464e953e47cSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
465e953e47cSHong Zhang       len         += len_si[proc];
466e953e47cSHong Zhang     }
467e953e47cSHong Zhang   }
468e953e47cSHong Zhang 
469e953e47cSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
470e953e47cSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
471e953e47cSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
472e953e47cSHong Zhang 
473e953e47cSHong Zhang   /* post the Irecv and Isend of coj */
474e953e47cSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
475e953e47cSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
476e953e47cSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
477e953e47cSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
478e953e47cSHong Zhang     if (!len_s[proc]) continue;
479e953e47cSHong Zhang     i    = owners_co[proc];
480e953e47cSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
481e953e47cSHong Zhang     k++;
482e953e47cSHong Zhang   }
483e953e47cSHong Zhang 
484e953e47cSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
485e953e47cSHong Zhang   /* ---------------------------------------- */
486e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
487e953e47cSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
488e953e47cSHong Zhang 
489e953e47cSHong Zhang   /* receives coj are complete */
490e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
491e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
492e953e47cSHong Zhang   }
493e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
494e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
495e953e47cSHong Zhang 
496e953e47cSHong Zhang   /* add received column indices into ta to update Crmax */
497e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
498e953e47cSHong Zhang     Jptr = buf_rj[k];
499e953e47cSHong Zhang     for (j=0; j<len_r[k]; j++) {
500e953e47cSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
501e953e47cSHong Zhang     }
502e953e47cSHong Zhang   }
503e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
504e953e47cSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
505e953e47cSHong Zhang 
506e953e47cSHong Zhang   /* (4) send and recv coi */
507e953e47cSHong Zhang   /*-----------------------*/
508e953e47cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
509e953e47cSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
510e953e47cSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
511e953e47cSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
512e953e47cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
513e953e47cSHong Zhang     if (!len_s[proc]) continue;
514e953e47cSHong Zhang     /* form outgoing message for i-structure:
515e953e47cSHong Zhang          buf_si[0]:                 nrows to be sent
516e953e47cSHong Zhang                [1:nrows]:           row index (global)
517e953e47cSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
518e953e47cSHong Zhang     */
519e953e47cSHong Zhang     /*-------------------------------------------*/
520e953e47cSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
521e953e47cSHong Zhang     buf_si_i    = buf_si + nrows+1;
522e953e47cSHong Zhang     buf_si[0]   = nrows;
523e953e47cSHong Zhang     buf_si_i[0] = 0;
524e953e47cSHong Zhang     nrows       = 0;
525e953e47cSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
526e953e47cSHong Zhang       nzi = coi[i+1] - coi[i];
527e953e47cSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
528e953e47cSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
529e953e47cSHong Zhang       nrows++;
530e953e47cSHong Zhang     }
531e953e47cSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
532e953e47cSHong Zhang     k++;
533e953e47cSHong Zhang     buf_si += len_si[proc];
534e953e47cSHong Zhang   }
535e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
536e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
537e953e47cSHong Zhang   }
538e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
539e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
540e953e47cSHong Zhang 
541e953e47cSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
542e953e47cSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
543e953e47cSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
544e953e47cSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
545b4e8d1b6SHong Zhang 
546e953e47cSHong Zhang   /* (5) compute the local portion of Cmpi      */
547e953e47cSHong Zhang   /* ------------------------------------------ */
548e953e47cSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
549e953e47cSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
550e953e47cSHong Zhang   current_space = free_space;
551e953e47cSHong Zhang 
552e953e47cSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
553e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {
554e953e47cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
555e953e47cSHong Zhang     nrows       = *buf_ri_k[k];
556e953e47cSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
557e953e47cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
558e953e47cSHong Zhang   }
559e953e47cSHong Zhang 
560e953e47cSHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
561cf742fccSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
562e953e47cSHong Zhang   for (i=0; i<pn; i++) {
563e953e47cSHong Zhang     /* add C_loc into Cmpi */
564e953e47cSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
565e953e47cSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
566cf742fccSHong Zhang     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
567e953e47cSHong Zhang 
568e953e47cSHong Zhang     /* add received col data into lnk */
569e953e47cSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
570e953e47cSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
571e953e47cSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
572e953e47cSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
573cf742fccSHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
574e953e47cSHong Zhang         nextrow[k]++; nextci[k]++;
575e953e47cSHong Zhang       }
576e953e47cSHong Zhang     }
577e953e47cSHong Zhang     nzi = lnk[0];
578e953e47cSHong Zhang 
579e953e47cSHong Zhang     /* copy data into free space, then initialize lnk */
580cf742fccSHong Zhang     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
581e953e47cSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
582e953e47cSHong Zhang   }
583e953e47cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
584cf742fccSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
585e953e47cSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
586e953e47cSHong Zhang 
587e953e47cSHong Zhang   /* local sizes and preallocation */
588e953e47cSHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
589e953e47cSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
590e953e47cSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
591e953e47cSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
592e953e47cSHong Zhang 
593e953e47cSHong Zhang   /* members in merge */
594e953e47cSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
595e953e47cSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
596e953e47cSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
597e953e47cSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
598e953e47cSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
599e953e47cSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
600e953e47cSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
601e953e47cSHong Zhang 
602e953e47cSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
603e953e47cSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
604e953e47cSHong Zhang   c->ptap         = ptap;
605e953e47cSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
606e953e47cSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
607cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
608e953e47cSHong Zhang 
609e953e47cSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
610e953e47cSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
611a4ffb656SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
612e953e47cSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
613cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
614e953e47cSHong Zhang   *C                     = Cmpi;
615e953e47cSHong Zhang   PetscFunctionReturn(0);
616e953e47cSHong Zhang }
617e953e47cSHong Zhang 
618e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
619e953e47cSHong Zhang {
620e953e47cSHong Zhang   PetscErrorCode      ierr;
621e953e47cSHong Zhang   Mat_PtAPMPI         *ptap;
622e953e47cSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
623e953e47cSHong Zhang   MPI_Comm            comm;
624e953e47cSHong Zhang   PetscMPIInt         size,rank;
625e953e47cSHong Zhang   Mat                 Cmpi;
626e953e47cSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
627e953e47cSHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
628e953e47cSHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
629e953e47cSHong Zhang   PetscBT             lnkbt;
630e953e47cSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
631e953e47cSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
632e953e47cSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
633e953e47cSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
634e953e47cSHong Zhang   MPI_Request         *swaits,*rwaits;
635e953e47cSHong Zhang   MPI_Status          *sstatus,rstatus;
636e953e47cSHong Zhang   PetscLayout         rowmap;
637e953e47cSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
638e953e47cSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
639e953e47cSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
640ec07b8f8SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
641e953e47cSHong Zhang   PetscScalar         *apv;
642e953e47cSHong Zhang   PetscTable          ta;
643b92f168fSBarry Smith   MatType             mtype;
644e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
645e953e47cSHong Zhang   PetscReal           apfill;
646e953e47cSHong Zhang #endif
647e953e47cSHong Zhang 
648e953e47cSHong Zhang   PetscFunctionBegin;
649b4e8d1b6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
650e953e47cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
651e953e47cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
652e953e47cSHong Zhang 
65352f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
654ec07b8f8SHong Zhang 
655e953e47cSHong Zhang   /* create symbolic parallel matrix Cmpi */
656e953e47cSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
657b92f168fSBarry Smith   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
658b92f168fSBarry Smith   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
6593697aca6SHong Zhang 
6603697aca6SHong Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-matptap_reuse",&Cmpi->reuse,NULL);CHKERRQ(ierr);
661e953e47cSHong Zhang 
662e953e47cSHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
663e953e47cSHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
664e953e47cSHong Zhang 
66515a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
66615a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
66715a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
668a4ffb656SHong Zhang   ptap->algType = 1;
66915a3b8e2SHong Zhang 
67015a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
67115a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
67215a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
67315a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
67415a3b8e2SHong Zhang 
67567a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
67667a12041SHong Zhang   /* --------------------------------- */
677419ecdd9Sandi selinger   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
678419ecdd9Sandi selinger   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
67915a3b8e2SHong Zhang 
68067a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
68167a12041SHong Zhang   /* ----------------------------------------------------------------- */
68267a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
68352f7967eSHong Zhang   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
684445158ffSHong Zhang 
6859a6dcab7SHong Zhang   /* create and initialize a linked list */
68645d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
6874b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
6884b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
68978d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
690d0e9a2f7SHong 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); */
69178d30b94SHong Zhang 
69278d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
693445158ffSHong Zhang 
6948cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
695ec07b8f8SHong Zhang   if (ao) {
696f91af8c7SBarry Smith     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
697ec07b8f8SHong Zhang   } else {
698ec07b8f8SHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
699ec07b8f8SHong Zhang   }
700445158ffSHong Zhang   current_space = free_space;
70167a12041SHong Zhang   nspacedouble  = 0;
70267a12041SHong Zhang 
70367a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
70467a12041SHong Zhang   api[0] = 0;
705445158ffSHong Zhang   for (i=0; i<am; i++) {
70667a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
70767a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
70867a12041SHong Zhang     nzi = ai[i+1] - ai[i];
70967a12041SHong Zhang     aj  = ad->j + ai[i];
710445158ffSHong Zhang     for (j=0; j<nzi; j++) {
711445158ffSHong Zhang       row  = aj[j];
71267a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
71367a12041SHong Zhang       Jptr = p_loc->j + pi[row];
714445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
715445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
716445158ffSHong Zhang     }
71767a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
718ec07b8f8SHong Zhang     if (ao) {
71967a12041SHong Zhang       ai = ao->i; pi = p_oth->i;
72067a12041SHong Zhang       nzi = ai[i+1] - ai[i];
72167a12041SHong Zhang       aj  = ao->j + ai[i];
722445158ffSHong Zhang       for (j=0; j<nzi; j++) {
723445158ffSHong Zhang         row  = aj[j];
72467a12041SHong Zhang         pnz  = pi[row+1] - pi[row];
72567a12041SHong Zhang         Jptr = p_oth->j + pi[row];
726445158ffSHong Zhang         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
727445158ffSHong Zhang       }
728ec07b8f8SHong Zhang     }
729445158ffSHong Zhang     apnz     = lnk[0];
730445158ffSHong Zhang     api[i+1] = api[i] + apnz;
731445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
732445158ffSHong Zhang 
733445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
734445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
735f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
736445158ffSHong Zhang       nspacedouble++;
737445158ffSHong Zhang     }
738445158ffSHong Zhang 
739445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
740445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
741445158ffSHong Zhang 
742445158ffSHong Zhang     current_space->array           += apnz;
743445158ffSHong Zhang     current_space->local_used      += apnz;
744445158ffSHong Zhang     current_space->local_remaining -= apnz;
745445158ffSHong Zhang   }
746681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
747445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
748445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
749445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
7509a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
7519a6dcab7SHong Zhang 
752aa690a28SHong Zhang   /* Create AP_loc for reuse */
753445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
754aa690a28SHong Zhang 
755aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
756ec07b8f8SHong Zhang   if (ao) {
757aa690a28SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
758ec07b8f8SHong Zhang   } else {
759ec07b8f8SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
760ec07b8f8SHong Zhang   }
761aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
762aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
763aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
764aa690a28SHong Zhang 
765aa690a28SHong Zhang   if (api[am]) {
766b11c1ec8SHong 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);
767aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
768aa690a28SHong Zhang   } else {
769b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
770aa690a28SHong Zhang   }
771aa690a28SHong Zhang #endif
772aa690a28SHong Zhang 
773681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
774681d504bSHong Zhang   /* ------------------------------------ */
77567a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
77615a3b8e2SHong Zhang 
77767a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
77867a12041SHong Zhang   /* ------------------------------------------ */
77915a3b8e2SHong Zhang   /* determine row ownership */
780445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
781445158ffSHong Zhang   rowmap->n  = pn;
782445158ffSHong Zhang   rowmap->bs = 1;
783445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
784445158ffSHong Zhang   owners = rowmap->range;
78515a3b8e2SHong Zhang 
78615a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
7878cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
7888cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
78915a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
79015a3b8e2SHong Zhang 
79167a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
79267a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
79367a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
79415a3b8e2SHong Zhang   proc  = 0;
79567a12041SHong Zhang   for (i=0; i<con; i++) {
79615a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
79715a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
79815a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
79915a3b8e2SHong Zhang   }
80015a3b8e2SHong Zhang 
80115a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
80215a3b8e2SHong Zhang   owners_co[0] = 0;
80367a12041SHong Zhang   nsend        = 0;
80415a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
80515a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
80615a3b8e2SHong Zhang     if (len_s[proc]) {
807445158ffSHong Zhang       nsend++;
80815a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
80915a3b8e2SHong Zhang       len         += len_si[proc];
81015a3b8e2SHong Zhang     }
81115a3b8e2SHong Zhang   }
81215a3b8e2SHong Zhang 
81315a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
814445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
815445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
81615a3b8e2SHong Zhang 
81715a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
81815a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
819445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
820445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
82115a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
82215a3b8e2SHong Zhang     if (!len_s[proc]) continue;
82315a3b8e2SHong Zhang     i    = owners_co[proc];
82415a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
82515a3b8e2SHong Zhang     k++;
82615a3b8e2SHong Zhang   }
82715a3b8e2SHong Zhang 
828681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
829681d504bSHong Zhang   /* ---------------------------------------- */
830681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
831681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
832681d504bSHong Zhang 
833681d504bSHong Zhang   /* receives coj are complete */
834445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
835445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
83615a3b8e2SHong Zhang   }
83715a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
838445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
83915a3b8e2SHong Zhang 
84078d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
84178d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
84278d30b94SHong Zhang     Jptr = buf_rj[k];
84378d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
84478d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
84578d30b94SHong Zhang     }
84678d30b94SHong Zhang   }
84778d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
84878d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
8499a6dcab7SHong Zhang 
85015a3b8e2SHong Zhang   /* (4) send and recv coi */
85115a3b8e2SHong Zhang   /*-----------------------*/
85215a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
853445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
85415a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
85515a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
85615a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
85715a3b8e2SHong Zhang     if (!len_s[proc]) continue;
85815a3b8e2SHong Zhang     /* form outgoing message for i-structure:
85915a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
86015a3b8e2SHong Zhang                [1:nrows]:           row index (global)
86115a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
86215a3b8e2SHong Zhang     */
86315a3b8e2SHong Zhang     /*-------------------------------------------*/
86415a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
86515a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
86615a3b8e2SHong Zhang     buf_si[0]   = nrows;
86715a3b8e2SHong Zhang     buf_si_i[0] = 0;
86815a3b8e2SHong Zhang     nrows       = 0;
86915a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
87015a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
87115a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
87215a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
87315a3b8e2SHong Zhang       nrows++;
87415a3b8e2SHong Zhang     }
87515a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
87615a3b8e2SHong Zhang     k++;
87715a3b8e2SHong Zhang     buf_si += len_si[proc];
87815a3b8e2SHong Zhang   }
879681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
880445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
88115a3b8e2SHong Zhang   }
88215a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
883445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
88415a3b8e2SHong Zhang 
8858cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
88615a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
88715a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
88815a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
889a4ffb656SHong Zhang 
89067a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
89167a12041SHong Zhang   /* ------------------------------------------ */
89278d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
89378d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
89415a3b8e2SHong Zhang   current_space = free_space;
89515a3b8e2SHong Zhang 
896445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
897445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
89815a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
89915a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
90015a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
90115a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
90215a3b8e2SHong Zhang   }
903445158ffSHong Zhang 
90478d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
90578d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
90615a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
90767a12041SHong Zhang     /* add C_loc into Cmpi */
90867a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
90967a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
91067a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
91115a3b8e2SHong Zhang 
91215a3b8e2SHong Zhang     /* add received col data into lnk */
913445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
91415a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
91515a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
91615a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
91715a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
91815a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
91915a3b8e2SHong Zhang       }
92015a3b8e2SHong Zhang     }
921d0e9a2f7SHong Zhang     nzi = lnk[0];
9228cb82516SHong Zhang 
92315a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
924d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
925d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
92615a3b8e2SHong Zhang   }
92715a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
92815a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
929445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
93080bb4639SHong Zhang 
931ae5f0867Sstefano_zampini   /* local sizes and preallocation */
93215a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
93315a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
93415a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
93515a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
93615a3b8e2SHong Zhang 
937445158ffSHong Zhang   /* members in merge */
938445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
939445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
940445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
941445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
942445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
943445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
944445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
94515a3b8e2SHong Zhang 
94615a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
94715a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
94815a3b8e2SHong Zhang   c->ptap         = ptap;
9491a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
9501a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
951cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
9528cb82516SHong Zhang   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
9532259aa2eSHong Zhang 
9541a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9551a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9561a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
957cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
9581a47ec54SHong Zhang   *C                     = Cmpi;
9590d3441aeSHong Zhang   PetscFunctionReturn(0);
9600d3441aeSHong Zhang }
9610d3441aeSHong Zhang 
962aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
9630d3441aeSHong Zhang {
9640d3441aeSHong Zhang   PetscErrorCode    ierr;
9650d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
9660d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
9672dd9e643SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
9680d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
9699ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
9700d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
9718cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
9728cb82516SHong Zhang   PetscScalar       *apa;
9730d3441aeSHong Zhang   const PetscInt    *cols;
9740d3441aeSHong Zhang   const PetscScalar *vals;
9750d3441aeSHong Zhang 
9760d3441aeSHong Zhang   PetscFunctionBegin;
977*a9899c97SHong Zhang   if (!C->reuse && ptap->reuse == MAT_REUSE_MATRIX) {
978*a9899c97SHong Zhang     MPI_Comm comm;
979*a9899c97SHong Zhang     ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
980*a9899c97SHong Zhang     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Call MatSetOption(C,MAT_REUSE,PETSC_TRUE) or '-matptap_reuse'");
981*a9899c97SHong Zhang   }
982*a9899c97SHong Zhang 
9830d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
984e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
985748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
986419ecdd9Sandi selinger     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
987419ecdd9Sandi selinger     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
988748c7196SHong Zhang   }
9890d3441aeSHong Zhang 
990e497d3c8SHong Zhang   /* 2) get AP_loc */
9910d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
9928cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
9930d3441aeSHong Zhang 
994e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
9950d3441aeSHong Zhang   /*-----------------------------------------------------*/
996748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
997748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
9980d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
9990d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1000e497d3c8SHong Zhang   }
10010d3441aeSHong Zhang 
10028cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
10038cb82516SHong Zhang   /* ---------------------------------------------- */
10040d3441aeSHong Zhang   /* get data from symbolic products */
10058cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
10062dd9e643SHong Zhang   if (ptap->P_oth) {
10078cb82516SHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
10082dd9e643SHong Zhang   }
10098cb82516SHong Zhang   apa   = ptap->apa;
1010681d504bSHong Zhang   api   = ap->i;
1011681d504bSHong Zhang   apj   = ap->j;
1012e497d3c8SHong Zhang   for (i=0; i<am; i++) {
1013445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
1014e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
1015e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
1016e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
1017e497d3c8SHong Zhang       col = apj[j+api[i]];
1018e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
1019e497d3c8SHong Zhang       apa[col] = 0.0;
10200d3441aeSHong Zhang     }
1021aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
10220d3441aeSHong Zhang   }
10230d3441aeSHong Zhang 
10248cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1025445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1026445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
10279ce11a7cSHong Zhang   C_loc = ptap->C_loc;
10289ce11a7cSHong Zhang   C_oth = ptap->C_oth;
10290d3441aeSHong Zhang 
10300d3441aeSHong Zhang   /* add C_loc and Co to to C */
10310d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
10320d3441aeSHong Zhang 
10330d3441aeSHong Zhang   /* C_loc -> C */
10340d3441aeSHong Zhang   cm    = C_loc->rmap->N;
10359ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
10368cb82516SHong Zhang   cols = c_seq->j;
10378cb82516SHong Zhang   vals = c_seq->a;
1038904d1e70Sandi selinger 
1039904d1e70Sandi selinger 
1040e9ede7d0Sandi selinger   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
1041904d1e70Sandi selinger   /* when there are no off-processor parts.  */
10421de21080Sandi selinger   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
10431de21080Sandi selinger   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
10441de21080Sandi selinger   /* a table, and other, more complex stuff has to be done. */
1045904d1e70Sandi selinger   if (C->assembled) {
1046904d1e70Sandi selinger     C->was_assembled = PETSC_TRUE;
1047904d1e70Sandi selinger     C->assembled     = PETSC_FALSE;
1048904d1e70Sandi selinger   }
1049904d1e70Sandi selinger   if (C->was_assembled) {
10500d3441aeSHong Zhang     for (i=0; i<cm; i++) {
10519ce11a7cSHong Zhang       ncols = c_seq->i[i+1] - c_seq->i[i];
10520d3441aeSHong Zhang       row = rstart + i;
1053904d1e70Sandi selinger       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10548cb82516SHong Zhang       cols += ncols; vals += ncols;
10550d3441aeSHong Zhang     }
1056904d1e70Sandi selinger   } else {
1057e9ede7d0Sandi selinger     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
1058904d1e70Sandi selinger   }
10590d3441aeSHong Zhang 
10600d3441aeSHong Zhang   /* Co -> C, off-processor part */
10619ce11a7cSHong Zhang   cm = C_oth->rmap->N;
10629ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
10638cb82516SHong Zhang   cols = c_seq->j;
10648cb82516SHong Zhang   vals = c_seq->a;
10659ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
10669ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
10670d3441aeSHong Zhang     row = p->garray[i];
10680d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10698cb82516SHong Zhang     cols += ncols; vals += ncols;
10700d3441aeSHong Zhang   }
1071904d1e70Sandi selinger 
10720d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10730d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10740d3441aeSHong Zhang 
1075748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
10763697aca6SHong Zhang 
10773697aca6SHong Zhang   if (!C->reuse) {
10783697aca6SHong Zhang     /* Free supporting data structure ptap */
10793697aca6SHong Zhang     ierr = MatDestroy_MPIAIJ_PtAPMPI_private(ptap);CHKERRQ(ierr);
10803697aca6SHong Zhang   }
10810d3441aeSHong Zhang   PetscFunctionReturn(0);
10820d3441aeSHong Zhang }
1083