xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision e83e5f86b8bcf4d7d2f46890b4075f5a0e7d07e4)
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;
193cdca5ebSHong Zhang   Mat_APMPI         *ptap=a->ap;
20cf97cf3cSHong Zhang   PetscBool         iascii;
21cf97cf3cSHong Zhang   PetscViewerFormat format;
22cf97cf3cSHong Zhang 
23cf97cf3cSHong Zhang   PetscFunctionBegin;
24baa1ba78SHong Zhang   if (!ptap) {
25baa1ba78SHong Zhang     /* hack: MatDuplicate() sets oldmat->ops->view to newmat which is a base mat class with null ptpa! */
26baa1ba78SHong Zhang     A->ops->view = MatView_MPIAIJ;
27baa1ba78SHong Zhang     ierr = (A->ops->view)(A,viewer);CHKERRQ(ierr);
28baa1ba78SHong Zhang     PetscFunctionReturn(0);
29baa1ba78SHong Zhang   }
30baa1ba78SHong Zhang 
31cf97cf3cSHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
32cf97cf3cSHong Zhang   if (iascii) {
33cf97cf3cSHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
34cf97cf3cSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
35cf97cf3cSHong Zhang       if (ptap->algType == 0) {
36cf97cf3cSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr);
37cf97cf3cSHong Zhang       } else if (ptap->algType == 1) {
38cf97cf3cSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr);
39cf97cf3cSHong Zhang       }
40cf97cf3cSHong Zhang     }
41cf97cf3cSHong Zhang   }
42cf97cf3cSHong Zhang   ierr = (ptap->view)(A,viewer);CHKERRQ(ierr);
43cf97cf3cSHong Zhang   PetscFunctionReturn(0);
44cf97cf3cSHong Zhang }
45cf97cf3cSHong Zhang 
464624976aSHong Zhang PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_AP(Mat A)
47cc6cb787SHong Zhang {
48cc6cb787SHong Zhang   PetscErrorCode      ierr;
49f8487c73SHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data;
503cdca5ebSHong Zhang   Mat_APMPI           *ptap=a->ap;
5160552ceaSHong Zhang   Mat_Merge_SeqsToMPI *merge;
52cc6cb787SHong Zhang 
53cc6cb787SHong Zhang   PetscFunctionBegin;
54dd4011a9SHong Zhang   if (!ptap) PetscFunctionReturn(0);
55dd4011a9SHong Zhang 
56b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
57f8487c73SHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
58a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
59a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
60c5af039cSHong Zhang   ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
6105d62848SHong Zhang   ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
6205d62848SHong Zhang   ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
638403a639SHong Zhang   if (ptap->AP_loc) { /* used by alg_rap */
64681d504bSHong Zhang     Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
65681d504bSHong Zhang     ierr = PetscFree(ap->i);CHKERRQ(ierr);
66681d504bSHong Zhang     ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
670d3441aeSHong Zhang     ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
688403a639SHong Zhang   } else { /* used by alg_ptap */
698403a639SHong Zhang     ierr = PetscFree(ptap->api);CHKERRQ(ierr);
708403a639SHong Zhang     ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
71681d504bSHong Zhang   }
722259aa2eSHong Zhang   ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
732259aa2eSHong Zhang   ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
74414894bdSHong Zhang   if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
75a560ef98SHong Zhang 
7660552ceaSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
7760552ceaSHong Zhang 
7860552ceaSHong Zhang   merge=ptap->merge;
798403a639SHong Zhang   if (merge) { /* used by alg_ptap */
80cc6cb787SHong Zhang     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
81cc6cb787SHong Zhang     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
82cc6cb787SHong Zhang     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
83cc6cb787SHong Zhang     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
84cc6cb787SHong Zhang     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
85c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
86cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
87c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
88cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
89445158ffSHong Zhang     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
90445158ffSHong Zhang     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
9105b42c5fSBarry Smith     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
926bf464f9SBarry Smith     ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
93f8487c73SHong Zhang     ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
94bf0cc555SLisandro Dalcin   }
95a3bb6f32SFande Kong   ierr = ISLocalToGlobalMappingDestroy(&ptap->ltog);CHKERRQ(ierr);
96a560ef98SHong Zhang 
97dd4011a9SHong Zhang   A->ops->destroy = ptap->destroy;
983cdca5ebSHong Zhang   ierr = PetscFree(a->ap);CHKERRQ(ierr);
993697aca6SHong Zhang   PetscFunctionReturn(0);
100c8b0795cSMark F. Adams }
101dce485f0SBarry Smith 
1023697aca6SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
1033697aca6SHong Zhang {
1043697aca6SHong Zhang   PetscErrorCode ierr;
1053697aca6SHong Zhang 
1063697aca6SHong Zhang   PetscFunctionBegin;
107dd4011a9SHong Zhang   ierr = (*A->ops->freeintermediatedatastructures)(A);CHKERRQ(ierr);
108dd4011a9SHong Zhang   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
109cc6cb787SHong Zhang   PetscFunctionReturn(0);
110cc6cb787SHong Zhang }
111cc6cb787SHong Zhang 
112150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
11342c57489SHong Zhang {
11442c57489SHong Zhang   PetscErrorCode ierr;
115a4ffb656SHong Zhang   PetscBool      flg;
11667a12041SHong Zhang   MPI_Comm       comm;
117a4ffb656SHong Zhang #if !defined(PETSC_HAVE_HYPRE)
118a4ffb656SHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
119a4ffb656SHong Zhang   PetscInt            nalg=2;
120a4ffb656SHong Zhang #else
121a4ffb656SHong Zhang   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
122a4ffb656SHong Zhang   PetscInt            nalg=3;
123a4ffb656SHong Zhang #endif
124a4ffb656SHong Zhang   PetscInt            pN=P->cmap->N,alg=1; /* set default algorithm */
12542c57489SHong Zhang 
12642c57489SHong Zhang   PetscFunctionBegin;
12767a12041SHong Zhang   /* check if matrix local sizes are compatible */
12867a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1296c4ed002SBarry 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);
1306c4ed002SBarry 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);
13167a12041SHong Zhang 
132cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
133a4ffb656SHong Zhang     /* pick an algorithm */
134715a5346SHong Zhang     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
135a4ffb656SHong Zhang     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
136a4ffb656SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
137a4ffb656SHong Zhang 
138a4ffb656SHong Zhang     if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */
139a4ffb656SHong Zhang       MatInfo     Ainfo,Pinfo;
140a4ffb656SHong Zhang       PetscInt    nz_local;
141a4ffb656SHong Zhang       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
142a4ffb656SHong Zhang 
143a4ffb656SHong Zhang       ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr);
144a4ffb656SHong Zhang       ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr);
145a4ffb656SHong Zhang       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
146a4ffb656SHong Zhang 
147a4ffb656SHong Zhang       if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
148a4ffb656SHong Zhang       ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr);
149a4ffb656SHong Zhang 
150a4ffb656SHong Zhang       if (alg_scalable) {
151a4ffb656SHong Zhang         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
1520d3441aeSHong Zhang       }
153a4ffb656SHong Zhang     }
154a4ffb656SHong Zhang 
155a4ffb656SHong Zhang     switch (alg) {
156a4ffb656SHong Zhang     case 1:
15780da8df7SHong Zhang       /* do R=P^T locally, then C=R*A*P -- nonscalable */
158a4ffb656SHong Zhang       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
159a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
1603ff4c91cSHong Zhang       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
161a4ffb656SHong Zhang       break;
162a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE)
163a4ffb656SHong Zhang     case 2:
164a4ffb656SHong Zhang       /* Use boomerAMGBuildCoarseOperator */
1650abfe76eSFande Kong       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
166a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
1670abfe76eSFande Kong       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
168a4ffb656SHong Zhang       break;
169a4ffb656SHong Zhang #endif
170a4ffb656SHong Zhang     default:
171a4ffb656SHong Zhang       /* do R=P^T locally, then C=R*A*P */
172a4ffb656SHong Zhang       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
173a4ffb656SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
174a4ffb656SHong Zhang       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
175a4ffb656SHong Zhang       break;
176a4ffb656SHong Zhang     }
1777d0a89b7SHong Zhang 
1786a85650dSHong Zhang     if (alg == 0 || alg == 1) {
1797d0a89b7SHong Zhang       Mat_MPIAIJ *c  = (Mat_MPIAIJ*)(*C)->data;
1807d0a89b7SHong Zhang       Mat_APMPI  *ap = c->ap;
1817d0a89b7SHong Zhang       ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr);
1827d0a89b7SHong Zhang       ap->freestruct = PETSC_FALSE;
1837d0a89b7SHong Zhang       ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr);
1847d0a89b7SHong Zhang       ierr = PetscOptionsEnd();CHKERRQ(ierr);
1857a7894deSKris Buschelman     }
1867d0a89b7SHong Zhang   }
1877d0a89b7SHong Zhang 
1883ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1898403a639SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
1903ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
19142c57489SHong Zhang   PetscFunctionReturn(0);
19242c57489SHong Zhang }
19342c57489SHong Zhang 
194cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C)
195cf742fccSHong Zhang {
196cf742fccSHong Zhang   PetscErrorCode    ierr;
197cf742fccSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
198cf742fccSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
199cf742fccSHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
2003cdca5ebSHong Zhang   Mat_APMPI         *ptap = c->ap;
201cf742fccSHong Zhang   Mat               AP_loc,C_loc,C_oth;
202a3bb6f32SFande Kong   PetscInt          i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz,nout;
203cf742fccSHong Zhang   PetscScalar       *apa;
204cf742fccSHong Zhang   const PetscInt    *cols;
205cf742fccSHong Zhang   const PetscScalar *vals;
206cf742fccSHong Zhang 
207cf742fccSHong Zhang   PetscFunctionBegin;
20880da8df7SHong Zhang   if (!ptap) {
20980da8df7SHong Zhang     MPI_Comm comm;
21080da8df7SHong Zhang     ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
21180da8df7SHong Zhang     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
21280da8df7SHong Zhang   }
21380da8df7SHong Zhang 
214cf742fccSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
215cf742fccSHong Zhang 
216cf742fccSHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
217cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
218419ecdd9Sandi selinger     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
219419ecdd9Sandi selinger     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
220cf742fccSHong Zhang   }
221cf742fccSHong Zhang 
222cf742fccSHong Zhang   /* 2) get AP_loc */
223cf742fccSHong Zhang   AP_loc = ptap->AP_loc;
224cf742fccSHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
225cf742fccSHong Zhang 
226cf742fccSHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
227cf742fccSHong Zhang   /*-----------------------------------------------------*/
228cf742fccSHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
229cf742fccSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
230cf742fccSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
231cf742fccSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
232cf742fccSHong Zhang   }
233cf742fccSHong Zhang 
234cf742fccSHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
235cf742fccSHong Zhang   /* ---------------------------------------------- */
236cf742fccSHong Zhang   /* get data from symbolic products */
237cf742fccSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
23852f7967eSHong Zhang   if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
239c072d3e2SSatish Balay   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!");
24052f7967eSHong Zhang 
241cf742fccSHong Zhang   api   = ap->i;
242cf742fccSHong Zhang   apj   = ap->j;
243a3bb6f32SFande Kong   ierr = ISLocalToGlobalMappingApply(ptap->ltog,api[AP_loc->rmap->n],apj,apj);CHKERRQ(ierr);
244cf742fccSHong Zhang   for (i=0; i<am; i++) {
245cf742fccSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
246cf742fccSHong Zhang     apnz = api[i+1] - api[i];
247b4e8d1b6SHong Zhang     apa = ap->a + api[i];
248b4e8d1b6SHong Zhang     ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr);
249b4e8d1b6SHong Zhang     AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa);
250cf742fccSHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
251cf742fccSHong Zhang   }
252a3bb6f32SFande Kong   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,api[AP_loc->rmap->n],apj,&nout,apj);CHKERRQ(ierr);
253a3bb6f32SFande Kong   if (api[AP_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",api[AP_loc->rmap->n],nout);
254cf742fccSHong Zhang 
255cf742fccSHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
256154d0d78SFande Kong   /* Always use scalable version since we are in the MPI scalable version */
257cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
258cf742fccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
259cf742fccSHong Zhang 
260cf742fccSHong Zhang   C_loc = ptap->C_loc;
261cf742fccSHong Zhang   C_oth = ptap->C_oth;
262cf742fccSHong Zhang 
263cf742fccSHong Zhang   /* add C_loc and Co to to C */
264cf742fccSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
265cf742fccSHong Zhang 
266cf742fccSHong Zhang   /* C_loc -> C */
267cf742fccSHong Zhang   cm    = C_loc->rmap->N;
268cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
269cf742fccSHong Zhang   cols = c_seq->j;
270cf742fccSHong Zhang   vals = c_seq->a;
271a3bb6f32SFande Kong   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_loc->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
272edeac6deSandi selinger 
273e9ede7d0Sandi selinger   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
274edeac6deSandi selinger   /* when there are no off-processor parts.  */
2751de21080Sandi selinger   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
2761de21080Sandi selinger   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
2771de21080Sandi selinger   /* a table, and other, more complex stuff has to be done. */
278edeac6deSandi selinger   if (C->assembled) {
279edeac6deSandi selinger     C->was_assembled = PETSC_TRUE;
280edeac6deSandi selinger     C->assembled     = PETSC_FALSE;
281edeac6deSandi selinger   }
282edeac6deSandi selinger   if (C->was_assembled) {
283cf742fccSHong Zhang     for (i=0; i<cm; i++) {
284cf742fccSHong Zhang       ncols = c_seq->i[i+1] - c_seq->i[i];
285cf742fccSHong Zhang       row = rstart + i;
286edeac6deSandi selinger       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
287cf742fccSHong Zhang       cols += ncols; vals += ncols;
288cf742fccSHong Zhang     }
289edeac6deSandi selinger   } else {
290e9ede7d0Sandi selinger     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
291edeac6deSandi selinger   }
292a3bb6f32SFande Kong   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_loc->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
293a3bb6f32SFande Kong   if (c_seq->i[C_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_seq->i[C_loc->rmap->n],nout);
294cf742fccSHong Zhang 
295cf742fccSHong Zhang   /* Co -> C, off-processor part */
296cf742fccSHong Zhang   cm = C_oth->rmap->N;
297cf742fccSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
298cf742fccSHong Zhang   cols = c_seq->j;
299cf742fccSHong Zhang   vals = c_seq->a;
300a3bb6f32SFande Kong   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_oth->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr);
301cf742fccSHong Zhang   for (i=0; i<cm; i++) {
302cf742fccSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
303cf742fccSHong Zhang     row = p->garray[i];
304cf742fccSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
305cf742fccSHong Zhang     cols += ncols; vals += ncols;
306cf742fccSHong Zhang   }
307cf742fccSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
308cf742fccSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
309cf742fccSHong Zhang 
310cf742fccSHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
31180da8df7SHong Zhang 
312a3bb6f32SFande Kong   ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_oth->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr);
313a3bb6f32SFande Kong   if (c_seq->i[C_oth->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_seq->i[C_loc->rmap->n],nout);
314a3bb6f32SFande Kong 
31580da8df7SHong Zhang   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
3167d0a89b7SHong Zhang   if (ptap->freestruct) {
31780da8df7SHong Zhang     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
31880da8df7SHong Zhang   }
319cf742fccSHong Zhang   PetscFunctionReturn(0);
320cf742fccSHong Zhang }
321cf742fccSHong Zhang 
322e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
3230d3441aeSHong Zhang {
3240d3441aeSHong Zhang   PetscErrorCode      ierr;
3253cdca5ebSHong Zhang   Mat_APMPI           *ptap;
326681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
3272259aa2eSHong Zhang   MPI_Comm            comm;
3282259aa2eSHong Zhang   PetscMPIInt         size,rank;
32976db11f6SHong Zhang   Mat                 Cmpi,P_loc,P_oth;
33015a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
331d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
332d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
333f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
33415a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
335d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
33615a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
33715a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
33815a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
339445158ffSHong Zhang   PetscLayout         rowmap;
340445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
341445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
342a3bb6f32SFande Kong   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi,nout;
34352f7967eSHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
34467a12041SHong Zhang   PetscScalar         *apv;
34578d30b94SHong Zhang   PetscTable          ta;
346b92f168fSBarry Smith   MatType             mtype;
347*e83e5f86SFande Kong   const char          *prefix;
348aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
3498cb82516SHong Zhang   PetscReal           apfill;
350aa690a28SHong Zhang #endif
35167a12041SHong Zhang 
35267a12041SHong Zhang   PetscFunctionBegin;
35367a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
35467a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
35567a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
356ae5f0867Sstefano_zampini 
35752f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
35852f7967eSHong Zhang 
359ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
360ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
361b92f168fSBarry Smith   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
362b92f168fSBarry Smith   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
363ae5f0867Sstefano_zampini 
3643cdca5ebSHong Zhang   /* create struct Mat_APMPI and attached it to C later */
365e953e47cSHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
366e953e47cSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
367cf97cf3cSHong Zhang   ptap->algType = 0;
368e953e47cSHong Zhang 
369e953e47cSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
37076db11f6SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr);
371e953e47cSHong Zhang   /* get P_loc by taking all local rows of P */
37276db11f6SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr);
37376db11f6SHong Zhang 
37476db11f6SHong Zhang   ptap->P_loc = P_loc;
37576db11f6SHong Zhang   ptap->P_oth = P_oth;
376e953e47cSHong Zhang 
377e953e47cSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
378e953e47cSHong Zhang   /* --------------------------------- */
379419ecdd9Sandi selinger   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
380419ecdd9Sandi selinger   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
381e953e47cSHong Zhang 
382e953e47cSHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
383e953e47cSHong Zhang   /* ----------------------------------------------------------------- */
38476db11f6SHong Zhang   p_loc  = (Mat_SeqAIJ*)P_loc->data;
38552f7967eSHong Zhang   if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data;
386e953e47cSHong Zhang 
387e953e47cSHong Zhang   /* create and initialize a linked list */
388e953e47cSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
38976db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta);
39076db11f6SHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta);
391e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
392e953e47cSHong Zhang 
39376db11f6SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
394e953e47cSHong Zhang 
395e953e47cSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
39652f7967eSHong Zhang   if (ao) {
397e953e47cSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
39852f7967eSHong Zhang   } else {
39952f7967eSHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
40052f7967eSHong Zhang   }
401e953e47cSHong Zhang   current_space = free_space;
402e953e47cSHong Zhang   nspacedouble  = 0;
403e953e47cSHong Zhang 
404e953e47cSHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
405e953e47cSHong Zhang   api[0] = 0;
406e953e47cSHong Zhang   for (i=0; i<am; i++) {
407e953e47cSHong Zhang     /* diagonal portion: Ad[i,:]*P */
408e953e47cSHong Zhang     ai = ad->i; pi = p_loc->i;
409e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
410e953e47cSHong Zhang     aj  = ad->j + ai[i];
411e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
412e953e47cSHong Zhang       row  = aj[j];
413e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
414e953e47cSHong Zhang       Jptr = p_loc->j + pi[row];
415e953e47cSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
41676db11f6SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
417e953e47cSHong Zhang     }
418e953e47cSHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
41952f7967eSHong Zhang     if (ao) {
420e953e47cSHong Zhang       ai = ao->i; pi = p_oth->i;
421e953e47cSHong Zhang       nzi = ai[i+1] - ai[i];
422e953e47cSHong Zhang       aj  = ao->j + ai[i];
423e953e47cSHong Zhang       for (j=0; j<nzi; j++) {
424e953e47cSHong Zhang         row  = aj[j];
425e953e47cSHong Zhang         pnz  = pi[row+1] - pi[row];
426e953e47cSHong Zhang         Jptr = p_oth->j + pi[row];
42776db11f6SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
428e953e47cSHong Zhang       }
42952f7967eSHong Zhang     }
430e953e47cSHong Zhang     apnz     = lnk[0];
431e953e47cSHong Zhang     api[i+1] = api[i] + apnz;
432e953e47cSHong Zhang 
433e953e47cSHong Zhang     /* if free space is not available, double the total space in the list */
434e953e47cSHong Zhang     if (current_space->local_remaining<apnz) {
435e953e47cSHong Zhang       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
436e953e47cSHong Zhang       nspacedouble++;
437e953e47cSHong Zhang     }
438e953e47cSHong Zhang 
439e953e47cSHong Zhang     /* Copy data into free space, then initialize lnk */
44076db11f6SHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
441e953e47cSHong Zhang 
442e953e47cSHong Zhang     current_space->array           += apnz;
443e953e47cSHong Zhang     current_space->local_used      += apnz;
444e953e47cSHong Zhang     current_space->local_remaining -= apnz;
445e953e47cSHong Zhang   }
446e953e47cSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
447e953e47cSHong Zhang   /* destroy list of free space and other temporary array(s) */
448a3bb6f32SFande Kong   ierr = PetscCalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
449e953e47cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
45076db11f6SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
451e953e47cSHong Zhang 
452e953e47cSHong Zhang   /* Create AP_loc for reuse */
453e953e47cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
454a3bb6f32SFande Kong   ierr = MatSeqAIJCompactOutExtraColumns_SeqAIJ(ptap->AP_loc, &ptap->ltog);CHKERRQ(ierr);
455e953e47cSHong Zhang 
456e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
45752f7967eSHong Zhang   if (ao) {
458e953e47cSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
45952f7967eSHong Zhang   } else {
46052f7967eSHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
46152f7967eSHong Zhang   }
462e953e47cSHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
463e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
464e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
465e953e47cSHong Zhang 
466e953e47cSHong Zhang   if (api[am]) {
467b11c1ec8SHong 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);
468e953e47cSHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
469e953e47cSHong Zhang   } else {
470b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
471e953e47cSHong Zhang   }
472e953e47cSHong Zhang #endif
473e953e47cSHong Zhang 
474e953e47cSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
475e953e47cSHong Zhang   /* ------------------------------------ */
476*e83e5f86SFande Kong   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
477*e83e5f86SFande Kong   ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr);
478*e83e5f86SFande Kong   ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
479e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
480e953e47cSHong Zhang 
481e953e47cSHong Zhang   /* (3) send coj of C_oth to other processors  */
482e953e47cSHong Zhang   /* ------------------------------------------ */
483e953e47cSHong Zhang   /* determine row ownership */
484e953e47cSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
485e953e47cSHong Zhang   rowmap->n  = pn;
486e953e47cSHong Zhang   rowmap->bs = 1;
487e953e47cSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
488e953e47cSHong Zhang   owners = rowmap->range;
489e953e47cSHong Zhang 
490e953e47cSHong Zhang   /* determine the number of messages to send, their lengths */
491e953e47cSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
492e953e47cSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
493e953e47cSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
494e953e47cSHong Zhang 
495e953e47cSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
496e953e47cSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
497e953e47cSHong Zhang   con   = ptap->C_oth->rmap->n;
498e953e47cSHong Zhang   proc  = 0;
499a3bb6f32SFande Kong   ierr = ISLocalToGlobalMappingApply(ptap->ltog,coi[con],coj,coj);CHKERRQ(ierr);
500e953e47cSHong Zhang   for (i=0; i<con; i++) {
501e953e47cSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
502e953e47cSHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
503e953e47cSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
504e953e47cSHong Zhang   }
505e953e47cSHong Zhang 
506e953e47cSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
507e953e47cSHong Zhang   owners_co[0] = 0;
508e953e47cSHong Zhang   nsend        = 0;
509e953e47cSHong Zhang   for (proc=0; proc<size; proc++) {
510e953e47cSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
511e953e47cSHong Zhang     if (len_s[proc]) {
512e953e47cSHong Zhang       nsend++;
513e953e47cSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
514e953e47cSHong Zhang       len         += len_si[proc];
515e953e47cSHong Zhang     }
516e953e47cSHong Zhang   }
517e953e47cSHong Zhang 
518e953e47cSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
519e953e47cSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
520e953e47cSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
521e953e47cSHong Zhang 
522e953e47cSHong Zhang   /* post the Irecv and Isend of coj */
523e953e47cSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
524e953e47cSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
525e953e47cSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
526e953e47cSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
527e953e47cSHong Zhang     if (!len_s[proc]) continue;
528e953e47cSHong Zhang     i    = owners_co[proc];
529e953e47cSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
530e953e47cSHong Zhang     k++;
531e953e47cSHong Zhang   }
532e953e47cSHong Zhang 
533e953e47cSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
534e953e47cSHong Zhang   /* ---------------------------------------- */
535*e83e5f86SFande Kong   ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr);
536*e83e5f86SFande Kong   ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
537e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
538e953e47cSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
539a3bb6f32SFande Kong   ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,c_loc->j);CHKERRQ(ierr);
540e953e47cSHong Zhang 
541e953e47cSHong Zhang   /* receives coj are complete */
542e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
543e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
544e953e47cSHong Zhang   }
545e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
546e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
547e953e47cSHong Zhang 
548e953e47cSHong Zhang   /* add received column indices into ta to update Crmax */
549e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
550e953e47cSHong Zhang     Jptr = buf_rj[k];
551e953e47cSHong Zhang     for (j=0; j<len_r[k]; j++) {
552e953e47cSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
553e953e47cSHong Zhang     }
554e953e47cSHong Zhang   }
555e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
556e953e47cSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
557e953e47cSHong Zhang 
558e953e47cSHong Zhang   /* (4) send and recv coi */
559e953e47cSHong Zhang   /*-----------------------*/
560e953e47cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
561e953e47cSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
562e953e47cSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
563e953e47cSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
564e953e47cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
565e953e47cSHong Zhang     if (!len_s[proc]) continue;
566e953e47cSHong Zhang     /* form outgoing message for i-structure:
567e953e47cSHong Zhang          buf_si[0]:                 nrows to be sent
568e953e47cSHong Zhang                [1:nrows]:           row index (global)
569e953e47cSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
570e953e47cSHong Zhang     */
571e953e47cSHong Zhang     /*-------------------------------------------*/
572e953e47cSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
573e953e47cSHong Zhang     buf_si_i    = buf_si + nrows+1;
574e953e47cSHong Zhang     buf_si[0]   = nrows;
575e953e47cSHong Zhang     buf_si_i[0] = 0;
576e953e47cSHong Zhang     nrows       = 0;
577e953e47cSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
578e953e47cSHong Zhang       nzi = coi[i+1] - coi[i];
579e953e47cSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
580e953e47cSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
581e953e47cSHong Zhang       nrows++;
582e953e47cSHong Zhang     }
583e953e47cSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
584e953e47cSHong Zhang     k++;
585e953e47cSHong Zhang     buf_si += len_si[proc];
586e953e47cSHong Zhang   }
587e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
588e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
589e953e47cSHong Zhang   }
590e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
591e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
592e953e47cSHong Zhang 
593e953e47cSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
594e953e47cSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
595e953e47cSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
596e953e47cSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
597b4e8d1b6SHong Zhang 
598e953e47cSHong Zhang   /* (5) compute the local portion of Cmpi      */
599e953e47cSHong Zhang   /* ------------------------------------------ */
600e953e47cSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
601e953e47cSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
602e953e47cSHong Zhang   current_space = free_space;
603e953e47cSHong Zhang 
604e953e47cSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
605e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {
606e953e47cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
607e953e47cSHong Zhang     nrows       = *buf_ri_k[k];
608e953e47cSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
609e953e47cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
610e953e47cSHong Zhang   }
611e953e47cSHong Zhang 
612e953e47cSHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
613cf742fccSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
614e953e47cSHong Zhang   for (i=0; i<pn; i++) {
615e953e47cSHong Zhang     /* add C_loc into Cmpi */
616e953e47cSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
617e953e47cSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
618cf742fccSHong Zhang     ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
619e953e47cSHong Zhang 
620e953e47cSHong Zhang     /* add received col data into lnk */
621e953e47cSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
622e953e47cSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
623e953e47cSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
624e953e47cSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
625cf742fccSHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
626e953e47cSHong Zhang         nextrow[k]++; nextci[k]++;
627e953e47cSHong Zhang       }
628e953e47cSHong Zhang     }
629e953e47cSHong Zhang     nzi = lnk[0];
630e953e47cSHong Zhang 
631e953e47cSHong Zhang     /* copy data into free space, then initialize lnk */
632cf742fccSHong Zhang     ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr);
633e953e47cSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
634e953e47cSHong Zhang   }
635e953e47cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
636cf742fccSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
637e953e47cSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
638e953e47cSHong Zhang 
639e953e47cSHong Zhang   /* local sizes and preallocation */
640e953e47cSHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
641e953e47cSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
642e953e47cSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
643e953e47cSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
644e953e47cSHong Zhang 
645e953e47cSHong Zhang   /* members in merge */
646e953e47cSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
647e953e47cSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
648e953e47cSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
649e953e47cSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
650e953e47cSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
651e953e47cSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
652e953e47cSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
653e953e47cSHong Zhang 
654e953e47cSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
655e953e47cSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
6563cdca5ebSHong Zhang   c->ap           = ptap;
657e953e47cSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
658e953e47cSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
659cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
660e953e47cSHong Zhang 
661e953e47cSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
662e953e47cSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
663a4ffb656SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable;
664e953e47cSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
665cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
6664624976aSHong Zhang   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
667e953e47cSHong Zhang   *C                     = Cmpi;
668a3bb6f32SFande Kong 
669a3bb6f32SFande Kong    nout = 0;
670a3bb6f32SFande Kong    ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_oth->i[ptap->C_oth->rmap->n],c_oth->j,&nout,c_oth->j);CHKERRQ(ierr);
671a3bb6f32SFande Kong    if (c_oth->i[ptap->C_oth->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_oth->i[ptap->C_oth->rmap->n],nout);
672a3bb6f32SFande Kong    ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,&nout,c_loc->j);CHKERRQ(ierr);
673a3bb6f32SFande Kong    if (c_loc->i[ptap->C_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_loc->i[ptap->C_loc->rmap->n],nout);
674a3bb6f32SFande Kong 
675e953e47cSHong Zhang   PetscFunctionReturn(0);
676e953e47cSHong Zhang }
677e953e47cSHong Zhang 
678e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
679e953e47cSHong Zhang {
680e953e47cSHong Zhang   PetscErrorCode      ierr;
6813cdca5ebSHong Zhang   Mat_APMPI           *ptap;
682e953e47cSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
683e953e47cSHong Zhang   MPI_Comm            comm;
684e953e47cSHong Zhang   PetscMPIInt         size,rank;
685e953e47cSHong Zhang   Mat                 Cmpi;
686e953e47cSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
687e953e47cSHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
688e953e47cSHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
689e953e47cSHong Zhang   PetscBT             lnkbt;
690e953e47cSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
691e953e47cSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
692e953e47cSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
693e953e47cSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
694e953e47cSHong Zhang   MPI_Request         *swaits,*rwaits;
695e953e47cSHong Zhang   MPI_Status          *sstatus,rstatus;
696e953e47cSHong Zhang   PetscLayout         rowmap;
697e953e47cSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
698e953e47cSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
699e953e47cSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
700ec07b8f8SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth;
701e953e47cSHong Zhang   PetscScalar         *apv;
702e953e47cSHong Zhang   PetscTable          ta;
703b92f168fSBarry Smith   MatType             mtype;
704*e83e5f86SFande Kong   const char          *prefix;
705e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
706e953e47cSHong Zhang   PetscReal           apfill;
707e953e47cSHong Zhang #endif
708e953e47cSHong Zhang 
709e953e47cSHong Zhang   PetscFunctionBegin;
710b4e8d1b6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
711e953e47cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
712e953e47cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
713e953e47cSHong Zhang 
71452f7967eSHong Zhang   if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data;
715ec07b8f8SHong Zhang 
716e953e47cSHong Zhang   /* create symbolic parallel matrix Cmpi */
717e953e47cSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
718b92f168fSBarry Smith   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
719b92f168fSBarry Smith   ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr);
720e953e47cSHong Zhang 
721e953e47cSHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
722e953e47cSHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
723e953e47cSHong Zhang 
7243cdca5ebSHong Zhang   /* create struct Mat_APMPI and attached it to C later */
72515a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
72615a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
727a4ffb656SHong Zhang   ptap->algType = 1;
72815a3b8e2SHong Zhang 
72915a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
73015a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
73115a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
73215a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
73315a3b8e2SHong Zhang 
73467a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
73567a12041SHong Zhang   /* --------------------------------- */
736419ecdd9Sandi selinger   ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
737419ecdd9Sandi selinger   ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
73815a3b8e2SHong Zhang 
73967a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
74067a12041SHong Zhang   /* ----------------------------------------------------------------- */
74167a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
74252f7967eSHong Zhang   if (ptap->P_oth) p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
743445158ffSHong Zhang 
7449a6dcab7SHong Zhang   /* create and initialize a linked list */
74545d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
7464b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
7474b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
74878d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
749d0e9a2f7SHong 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); */
75078d30b94SHong Zhang 
75178d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
752445158ffSHong Zhang 
7538cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
754ec07b8f8SHong Zhang   if (ao) {
755f91af8c7SBarry Smith     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
756ec07b8f8SHong Zhang   } else {
757ec07b8f8SHong Zhang     ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr);
758ec07b8f8SHong Zhang   }
759445158ffSHong Zhang   current_space = free_space;
76067a12041SHong Zhang   nspacedouble  = 0;
76167a12041SHong Zhang 
76267a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
76367a12041SHong Zhang   api[0] = 0;
764445158ffSHong Zhang   for (i=0; i<am; i++) {
76567a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
76667a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
76767a12041SHong Zhang     nzi = ai[i+1] - ai[i];
76867a12041SHong Zhang     aj  = ad->j + ai[i];
769445158ffSHong Zhang     for (j=0; j<nzi; j++) {
770445158ffSHong Zhang       row  = aj[j];
77167a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
77267a12041SHong Zhang       Jptr = p_loc->j + pi[row];
773445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
774445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
775445158ffSHong Zhang     }
77667a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
777ec07b8f8SHong Zhang     if (ao) {
77867a12041SHong Zhang       ai = ao->i; pi = p_oth->i;
77967a12041SHong Zhang       nzi = ai[i+1] - ai[i];
78067a12041SHong Zhang       aj  = ao->j + ai[i];
781445158ffSHong Zhang       for (j=0; j<nzi; j++) {
782445158ffSHong Zhang         row  = aj[j];
78367a12041SHong Zhang         pnz  = pi[row+1] - pi[row];
78467a12041SHong Zhang         Jptr = p_oth->j + pi[row];
785445158ffSHong Zhang         ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
786445158ffSHong Zhang       }
787ec07b8f8SHong Zhang     }
788445158ffSHong Zhang     apnz     = lnk[0];
789445158ffSHong Zhang     api[i+1] = api[i] + apnz;
790445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
791445158ffSHong Zhang 
792445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
793445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
794f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
795445158ffSHong Zhang       nspacedouble++;
796445158ffSHong Zhang     }
797445158ffSHong Zhang 
798445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
799445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
800445158ffSHong Zhang 
801445158ffSHong Zhang     current_space->array           += apnz;
802445158ffSHong Zhang     current_space->local_used      += apnz;
803445158ffSHong Zhang     current_space->local_remaining -= apnz;
804445158ffSHong Zhang   }
805681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
806445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
807445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
808445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
8099a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
8109a6dcab7SHong Zhang 
811aa690a28SHong Zhang   /* Create AP_loc for reuse */
812445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
813aa690a28SHong Zhang 
814aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
815ec07b8f8SHong Zhang   if (ao) {
816aa690a28SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
817ec07b8f8SHong Zhang   } else {
818ec07b8f8SHong Zhang     apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1);
819ec07b8f8SHong Zhang   }
820aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
821aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
822aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
823aa690a28SHong Zhang 
824aa690a28SHong Zhang   if (api[am]) {
825b11c1ec8SHong 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);
826aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
827aa690a28SHong Zhang   } else {
828b11c1ec8SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr);
829aa690a28SHong Zhang   }
830aa690a28SHong Zhang #endif
831aa690a28SHong Zhang 
832681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
833681d504bSHong Zhang   /* ------------------------------------ */
834*e83e5f86SFande Kong   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
835*e83e5f86SFande Kong   ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr);
836*e83e5f86SFande Kong   ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr);
83767a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
83815a3b8e2SHong Zhang 
83967a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
84067a12041SHong Zhang   /* ------------------------------------------ */
84115a3b8e2SHong Zhang   /* determine row ownership */
842445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
843445158ffSHong Zhang   rowmap->n  = pn;
844445158ffSHong Zhang   rowmap->bs = 1;
845445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
846445158ffSHong Zhang   owners = rowmap->range;
84715a3b8e2SHong Zhang 
84815a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
8498cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
8508cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
85115a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
85215a3b8e2SHong Zhang 
85367a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
85467a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
85567a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
85615a3b8e2SHong Zhang   proc  = 0;
85767a12041SHong Zhang   for (i=0; i<con; i++) {
85815a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
85915a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
86015a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
86115a3b8e2SHong Zhang   }
86215a3b8e2SHong Zhang 
86315a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
86415a3b8e2SHong Zhang   owners_co[0] = 0;
86567a12041SHong Zhang   nsend        = 0;
86615a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
86715a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
86815a3b8e2SHong Zhang     if (len_s[proc]) {
869445158ffSHong Zhang       nsend++;
87015a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
87115a3b8e2SHong Zhang       len         += len_si[proc];
87215a3b8e2SHong Zhang     }
87315a3b8e2SHong Zhang   }
87415a3b8e2SHong Zhang 
87515a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
876445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
877445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
87815a3b8e2SHong Zhang 
87915a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
88015a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
881445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
882445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
88315a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
88415a3b8e2SHong Zhang     if (!len_s[proc]) continue;
88515a3b8e2SHong Zhang     i    = owners_co[proc];
88615a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
88715a3b8e2SHong Zhang     k++;
88815a3b8e2SHong Zhang   }
88915a3b8e2SHong Zhang 
890681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
891681d504bSHong Zhang   /* ---------------------------------------- */
892*e83e5f86SFande Kong   ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr);
893*e83e5f86SFande Kong   ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr);
894681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
895681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
896681d504bSHong Zhang 
897681d504bSHong Zhang   /* receives coj are complete */
898445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
899445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
90015a3b8e2SHong Zhang   }
90115a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
902445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
90315a3b8e2SHong Zhang 
90478d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
90578d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
90678d30b94SHong Zhang     Jptr = buf_rj[k];
90778d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
90878d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
90978d30b94SHong Zhang     }
91078d30b94SHong Zhang   }
91178d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
91278d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
9139a6dcab7SHong Zhang 
91415a3b8e2SHong Zhang   /* (4) send and recv coi */
91515a3b8e2SHong Zhang   /*-----------------------*/
91615a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
917445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
91815a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
91915a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
92015a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
92115a3b8e2SHong Zhang     if (!len_s[proc]) continue;
92215a3b8e2SHong Zhang     /* form outgoing message for i-structure:
92315a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
92415a3b8e2SHong Zhang                [1:nrows]:           row index (global)
92515a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
92615a3b8e2SHong Zhang     */
92715a3b8e2SHong Zhang     /*-------------------------------------------*/
92815a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
92915a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
93015a3b8e2SHong Zhang     buf_si[0]   = nrows;
93115a3b8e2SHong Zhang     buf_si_i[0] = 0;
93215a3b8e2SHong Zhang     nrows       = 0;
93315a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
93415a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
93515a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
93615a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
93715a3b8e2SHong Zhang       nrows++;
93815a3b8e2SHong Zhang     }
93915a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
94015a3b8e2SHong Zhang     k++;
94115a3b8e2SHong Zhang     buf_si += len_si[proc];
94215a3b8e2SHong Zhang   }
943681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
944445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
94515a3b8e2SHong Zhang   }
94615a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
947445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
94815a3b8e2SHong Zhang 
9498cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
95015a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
95115a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
95215a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
953a4ffb656SHong Zhang 
95467a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
95567a12041SHong Zhang   /* ------------------------------------------ */
95678d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
95778d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
95815a3b8e2SHong Zhang   current_space = free_space;
95915a3b8e2SHong Zhang 
960445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
961445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
96215a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
96315a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
96415a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
96515a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
96615a3b8e2SHong Zhang   }
967445158ffSHong Zhang 
96878d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
96978d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
97015a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
97167a12041SHong Zhang     /* add C_loc into Cmpi */
97267a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
97367a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
97467a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
97515a3b8e2SHong Zhang 
97615a3b8e2SHong Zhang     /* add received col data into lnk */
977445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
97815a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
97915a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
98015a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
98115a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
98215a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
98315a3b8e2SHong Zhang       }
98415a3b8e2SHong Zhang     }
985d0e9a2f7SHong Zhang     nzi = lnk[0];
9868cb82516SHong Zhang 
98715a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
988d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
989d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
99015a3b8e2SHong Zhang   }
99115a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
99215a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
993445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
99480bb4639SHong Zhang 
995ae5f0867Sstefano_zampini   /* local sizes and preallocation */
99615a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
99715a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
99815a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
99915a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
100015a3b8e2SHong Zhang 
1001445158ffSHong Zhang   /* members in merge */
1002445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
1003445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
1004445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
1005445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
1006445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
1007445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
1008445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
100915a3b8e2SHong Zhang 
101015a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
101115a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
10123cdca5ebSHong Zhang   c->ap           = ptap;
10131a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
10141a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1015cf97cf3cSHong Zhang   ptap->view      = Cmpi->ops->view;
10168cb82516SHong Zhang   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
10172259aa2eSHong Zhang 
10181a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
10191a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
10201a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1021cf97cf3cSHong Zhang   Cmpi->ops->view        = MatView_MPIAIJ_PtAP;
10224624976aSHong Zhang   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
10231a47ec54SHong Zhang   *C                     = Cmpi;
10240d3441aeSHong Zhang   PetscFunctionReturn(0);
10250d3441aeSHong Zhang }
10260d3441aeSHong Zhang 
1027aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
10280d3441aeSHong Zhang {
10290d3441aeSHong Zhang   PetscErrorCode    ierr;
10300d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
10310d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10322dd9e643SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth=NULL,*c_seq;
10333cdca5ebSHong Zhang   Mat_APMPI         *ptap = c->ap;
10349ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
10350d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
10368cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
10378cb82516SHong Zhang   PetscScalar       *apa;
10380d3441aeSHong Zhang   const PetscInt    *cols;
10390d3441aeSHong Zhang   const PetscScalar *vals;
10400d3441aeSHong Zhang 
10410d3441aeSHong Zhang   PetscFunctionBegin;
1042dd4011a9SHong Zhang   if (!ptap) {
1043a9899c97SHong Zhang     MPI_Comm comm;
1044a9899c97SHong Zhang     ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1045dd4011a9SHong Zhang     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1046a9899c97SHong Zhang   }
1047a9899c97SHong Zhang 
10480d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
1049e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
1050748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
1051419ecdd9Sandi selinger     ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
1052419ecdd9Sandi selinger     ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1053748c7196SHong Zhang   }
10540d3441aeSHong Zhang 
1055e497d3c8SHong Zhang   /* 2) get AP_loc */
10560d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
10578cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
10580d3441aeSHong Zhang 
1059e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
10600d3441aeSHong Zhang   /*-----------------------------------------------------*/
1061748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
1062748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
10630d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
10640d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1065e497d3c8SHong Zhang   }
10660d3441aeSHong Zhang 
10678cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
10688cb82516SHong Zhang   /* ---------------------------------------------- */
10690d3441aeSHong Zhang   /* get data from symbolic products */
10708cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
10712dd9e643SHong Zhang   if (ptap->P_oth) {
10728cb82516SHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
10732dd9e643SHong Zhang   }
10748cb82516SHong Zhang   apa   = ptap->apa;
1075681d504bSHong Zhang   api   = ap->i;
1076681d504bSHong Zhang   apj   = ap->j;
1077e497d3c8SHong Zhang   for (i=0; i<am; i++) {
1078445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
1079e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
1080e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
1081e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
1082e497d3c8SHong Zhang       col = apj[j+api[i]];
1083e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
1084e497d3c8SHong Zhang       apa[col] = 0.0;
10850d3441aeSHong Zhang     }
1086aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
10870d3441aeSHong Zhang   }
10880d3441aeSHong Zhang 
10898cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
1090445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
1091445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
10929ce11a7cSHong Zhang   C_loc = ptap->C_loc;
10939ce11a7cSHong Zhang   C_oth = ptap->C_oth;
10940d3441aeSHong Zhang 
10950d3441aeSHong Zhang   /* add C_loc and Co to to C */
10960d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
10970d3441aeSHong Zhang 
10980d3441aeSHong Zhang   /* C_loc -> C */
10990d3441aeSHong Zhang   cm    = C_loc->rmap->N;
11009ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
11018cb82516SHong Zhang   cols = c_seq->j;
11028cb82516SHong Zhang   vals = c_seq->a;
1103904d1e70Sandi selinger 
1104904d1e70Sandi selinger 
1105e9ede7d0Sandi selinger   /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */
1106904d1e70Sandi selinger   /* when there are no off-processor parts.  */
11071de21080Sandi selinger   /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */
11081de21080Sandi selinger   /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */
11091de21080Sandi selinger   /* a table, and other, more complex stuff has to be done. */
1110904d1e70Sandi selinger   if (C->assembled) {
1111904d1e70Sandi selinger     C->was_assembled = PETSC_TRUE;
1112904d1e70Sandi selinger     C->assembled     = PETSC_FALSE;
1113904d1e70Sandi selinger   }
1114904d1e70Sandi selinger   if (C->was_assembled) {
11150d3441aeSHong Zhang     for (i=0; i<cm; i++) {
11169ce11a7cSHong Zhang       ncols = c_seq->i[i+1] - c_seq->i[i];
11170d3441aeSHong Zhang       row = rstart + i;
1118904d1e70Sandi selinger       ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
11198cb82516SHong Zhang       cols += ncols; vals += ncols;
11200d3441aeSHong Zhang     }
1121904d1e70Sandi selinger   } else {
1122e9ede7d0Sandi selinger     ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr);
1123904d1e70Sandi selinger   }
11240d3441aeSHong Zhang 
11250d3441aeSHong Zhang   /* Co -> C, off-processor part */
11269ce11a7cSHong Zhang   cm = C_oth->rmap->N;
11279ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
11288cb82516SHong Zhang   cols = c_seq->j;
11298cb82516SHong Zhang   vals = c_seq->a;
11309ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
11319ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
11320d3441aeSHong Zhang     row = p->garray[i];
11330d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
11348cb82516SHong Zhang     cols += ncols; vals += ncols;
11350d3441aeSHong Zhang   }
1136904d1e70Sandi selinger 
11370d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11380d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11390d3441aeSHong Zhang 
1140748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
11413697aca6SHong Zhang 
1142dd4011a9SHong Zhang   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
11437d0a89b7SHong Zhang   if (ptap->freestruct) {
1144dd4011a9SHong Zhang     ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
11453697aca6SHong Zhang   }
11460d3441aeSHong Zhang   PetscFunctionReturn(0);
11470d3441aeSHong Zhang }
1148