xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 96b95a6bf9ed541503cb02ce9899069044d5ddfa)
1be1d678aSKris Buschelman 
22499ec78SHong Zhang /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
42499ec78SHong Zhang           C = A * B
52499ec78SHong Zhang */
6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
8c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>
111634b032SHong Zhang 
122499ec78SHong Zhang #undef __FUNCT__
132499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
142499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
152499ec78SHong Zhang {
162499ec78SHong Zhang   PetscErrorCode ierr;
172499ec78SHong Zhang 
182499ec78SHong Zhang   PetscFunctionBegin;
192499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
20598bc09dSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
21a1a4e84aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
22598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
232499ec78SHong Zhang   }
24598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
25598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
26598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
272499ec78SHong Zhang   PetscFunctionReturn(0);
282499ec78SHong Zhang }
292499ec78SHong Zhang 
30f33d1a9aSHong Zhang #undef __FUNCT__
31a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
32a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
33598bc09dSHong Zhang {
34598bc09dSHong Zhang   PetscErrorCode ierr;
35598bc09dSHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
36598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
37598bc09dSHong Zhang 
38598bc09dSHong Zhang   PetscFunctionBegin;
39b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
40598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
41a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
42a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
43a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
44a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
45598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
46598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
47598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
48598bc09dSHong Zhang   PetscFunctionReturn(0);
49598bc09dSHong Zhang }
50598bc09dSHong Zhang 
512499ec78SHong Zhang #undef __FUNCT__
52a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
53a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
544ae0847bSHong Zhang {
554ae0847bSHong Zhang   PetscErrorCode     ierr;
564ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
574ae0847bSHong Zhang   Mat_PtAPMPI        *ptap=a->ptap;
584ae0847bSHong Zhang 
594ae0847bSHong Zhang   PetscFunctionBegin;
604ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
614ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
624ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
634ae0847bSHong Zhang   PetscFunctionReturn(0);
644ae0847bSHong Zhang }
654ae0847bSHong Zhang 
664ae0847bSHong Zhang #undef __FUNCT__
67a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
68a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
69598bc09dSHong Zhang {
70598bc09dSHong Zhang   PetscErrorCode     ierr;
714ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
72598bc09dSHong Zhang   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
73598bc09dSHong Zhang   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
74598bc09dSHong Zhang   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
75598bc09dSHong Zhang   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
76598bc09dSHong Zhang   Mat_SeqAIJ         *p_loc,*p_oth;
77598bc09dSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
78598bc09dSHong Zhang   PetscScalar        *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
79598bc09dSHong Zhang   PetscInt           cm=C->rmap->n,anz,pnz;
80598bc09dSHong Zhang   Mat_PtAPMPI        *ptap=c->ptap;
8165a57a32SSean Farley   PetscInt           *api,*apj,*apJ,i,j,k,row;
8229825010SHong Zhang   PetscInt           cstart=C->cmap->rstart;
83598bc09dSHong Zhang   PetscInt           cdnz,conz,k0,k1;
84598bc09dSHong Zhang 
85598bc09dSHong Zhang   PetscFunctionBegin;
86a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
87598bc09dSHong Zhang   /*-----------------------------------------------------*/
88a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
89b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
90a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
91598bc09dSHong Zhang 
92598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
93598bc09dSHong Zhang   /*----------------------------------------------------------*/
94598bc09dSHong Zhang   /* get data from symbolic products */
95a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
96a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
97598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
98598bc09dSHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
99598bc09dSHong Zhang 
100598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
101598bc09dSHong Zhang   apa = ptap->apa;
102598bc09dSHong Zhang 
10329825010SHong Zhang   api = ptap->api;
10429825010SHong Zhang   apj = ptap->apj;
105598bc09dSHong Zhang   for (i=0; i<cm; i++) {
106598bc09dSHong Zhang     /* diagonal portion of A */
107598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
108598bc09dSHong Zhang     adj = ad->j + adi[i];
109598bc09dSHong Zhang     ada = ad->a + adi[i];
110598bc09dSHong Zhang     for (j=0; j<anz; j++) {
111598bc09dSHong Zhang       row = adj[j];
112598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
113598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
114598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
115598bc09dSHong Zhang 
116598bc09dSHong Zhang       /* perform dense axpy */
117598bc09dSHong Zhang       valtmp = ada[j];
118598bc09dSHong Zhang       for (k=0; k<pnz; k++){
119598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
120598bc09dSHong Zhang       }
121598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
122598bc09dSHong Zhang     }
123598bc09dSHong Zhang 
124598bc09dSHong Zhang     /* off-diagonal portion of A */
125598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
126598bc09dSHong Zhang     aoj = ao->j + aoi[i];
127598bc09dSHong Zhang     aoa = ao->a + aoi[i];
128598bc09dSHong Zhang     for (j=0; j<anz; j++) {
129598bc09dSHong Zhang       row = aoj[j];
130598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
131598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
132598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
133598bc09dSHong Zhang 
134598bc09dSHong Zhang       /* perform dense axpy */
135598bc09dSHong Zhang       valtmp = aoa[j];
136598bc09dSHong Zhang       for (k=0; k<pnz; k++){
137598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
138598bc09dSHong Zhang       }
139598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
140598bc09dSHong Zhang     }
141598bc09dSHong Zhang 
142598bc09dSHong Zhang     /* set values in C */
143598bc09dSHong Zhang     apJ = apj + api[i];
144598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
145598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
146598bc09dSHong Zhang 
147598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
148598bc09dSHong Zhang     ca = coa + co->i[i];
149598bc09dSHong Zhang     k  = 0;
150598bc09dSHong Zhang     for (k0=0; k0<conz; k0++){
151598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
152598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
153598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
154598bc09dSHong Zhang       k++;
155598bc09dSHong Zhang     }
156598bc09dSHong Zhang 
157598bc09dSHong Zhang     /* diagonal part of C */
158598bc09dSHong Zhang     ca = cda + cd->i[i];
159598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++){
160598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
161598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
162598bc09dSHong Zhang       k++;
163598bc09dSHong Zhang     }
164598bc09dSHong Zhang 
165598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
166598bc09dSHong Zhang     ca = coa + co->i[i];
167598bc09dSHong Zhang     for (; k0<conz; k0++){
168598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
169598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
170598bc09dSHong Zhang       k++;
171598bc09dSHong Zhang     }
172598bc09dSHong Zhang   }
173598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
174598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175598bc09dSHong Zhang   PetscFunctionReturn(0);
176598bc09dSHong Zhang }
177598bc09dSHong Zhang 
178598bc09dSHong Zhang #undef __FUNCT__
179a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
180a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
181598bc09dSHong Zhang {
182598bc09dSHong Zhang   PetscErrorCode       ierr;
183598bc09dSHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
184bfcd1a73SHong Zhang   Mat                  Cmpi;
185598bc09dSHong Zhang   Mat_PtAPMPI          *ptap;
186598bc09dSHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
1874ae0847bSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
188bfcd1a73SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1894ae0847bSHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
1904ae0847bSHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
191d14fa243SHong Zhang   PetscInt             *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
192bfcd1a73SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
193598bc09dSHong Zhang   PetscBT              lnkbt;
194598bc09dSHong Zhang   PetscScalar          *apa;
195bfcd1a73SHong Zhang   PetscReal            afill;
196cf1a0672SHong Zhang   PetscBool            scalable=PETSC_FALSE;
197f84c536bSHong Zhang   PetscInt             nlnk_max,armax,prmax;
198598bc09dSHong Zhang 
199598bc09dSHong Zhang   PetscFunctionBegin;
200a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
201cf1a0672SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
202a1a4e84aSHong Zhang   }
203152983bfSHong Zhang 
204152983bfSHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
205cf1a0672SHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
206cf1a0672SHong Zhang     if (scalable){
20725023028SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr);
2081634b032SHong Zhang       PetscFunctionReturn(0);
2091634b032SHong Zhang     }
210152983bfSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
211a1a4e84aSHong Zhang 
212598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
213598bc09dSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
214598bc09dSHong Zhang 
215598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
216b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
21719f0e57aSHong Zhang 
218598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
219a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
220598bc09dSHong Zhang 
221a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
222a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
223598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
224598bc09dSHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
225598bc09dSHong Zhang 
226598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
227598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
228598bc09dSHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
229a1a4e84aSHong Zhang   ptap->api = api;
230598bc09dSHong Zhang   api[0]    = 0;
231598bc09dSHong Zhang 
232598bc09dSHong Zhang   /* create and initialize a linked list */
233f84c536bSHong Zhang   armax = ad->rmax+ao->rmax;
234f84c536bSHong Zhang   prmax = PetscMax(p_loc->rmax,p_oth->rmax);
235f84c536bSHong Zhang   nlnk_max = armax*prmax;
236f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2370d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
238598bc09dSHong Zhang 
239bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
240bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
241598bc09dSHong Zhang   current_space = free_space;
242598bc09dSHong Zhang 
243598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
244598bc09dSHong Zhang   for (i=0; i<am; i++) {
245598bc09dSHong Zhang     apnz = 0;
246598bc09dSHong Zhang     /* diagonal portion of A */
247598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
248598bc09dSHong Zhang     for (j=0; j<nzi; j++){
249598bc09dSHong Zhang       row = *adj++;
250598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
251598bc09dSHong Zhang       Jptr  = pj_loc + pi_loc[row];
252598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2531fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
254598bc09dSHong Zhang     }
255598bc09dSHong Zhang     /* off-diagonal portion of A */
256598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
257598bc09dSHong Zhang     for (j=0; j<nzi; j++){
258598bc09dSHong Zhang       row = *aoj++;
259598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
260598bc09dSHong Zhang       Jptr  = pj_oth + pi_oth[row];
2611fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
262598bc09dSHong Zhang     }
263598bc09dSHong Zhang 
264d14fa243SHong Zhang     apnz     = lnk[0];
265598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
266598bc09dSHong Zhang 
267598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
268598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
269598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
270598bc09dSHong Zhang       nspacedouble++;
271598bc09dSHong Zhang     }
272598bc09dSHong Zhang 
273598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2741fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
275598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
276598bc09dSHong Zhang     current_space->array           += apnz;
277598bc09dSHong Zhang     current_space->local_used      += apnz;
278598bc09dSHong Zhang     current_space->local_remaining -= apnz;
279598bc09dSHong Zhang   }
280598bc09dSHong Zhang 
281598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
282598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
283a1a4e84aSHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
284a1a4e84aSHong Zhang   apj  = ptap->apj;
285a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
286598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
287598bc09dSHong Zhang 
288f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
289f84c536bSHong Zhang   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
290f84c536bSHong Zhang   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
291f84c536bSHong Zhang   ptap->apa = apa;
292f84c536bSHong Zhang 
293bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
294598bc09dSHong Zhang   /*----------------------------------------------------*/
295bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
296bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
297bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
298bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
299598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
300598bc09dSHong Zhang   for (i=0; i<am; i++){
301598bc09dSHong Zhang     row  = i + rstart;
302598bc09dSHong Zhang     apnz = api[i+1] - api[i];
303bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
304598bc09dSHong Zhang     apj += apnz;
305598bc09dSHong Zhang   }
306bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
307bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
308598bc09dSHong Zhang 
309bfcd1a73SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
310bfcd1a73SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
311bfcd1a73SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
312bfcd1a73SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
313bfcd1a73SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
314598bc09dSHong Zhang 
315bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
316bfcd1a73SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
317598bc09dSHong Zhang   c->ptap  = ptap;
318598bc09dSHong Zhang 
319bfcd1a73SHong Zhang   *C = Cmpi;
320bfcd1a73SHong Zhang 
321bfcd1a73SHong Zhang   /* set MatInfo */
322bfcd1a73SHong Zhang   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
323bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
324bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
325bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
326bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
327bfcd1a73SHong Zhang 
328bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
329bfcd1a73SHong Zhang   if (api[am]) {
330bfcd1a73SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
331bfcd1a73SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
332bfcd1a73SHong Zhang   } else {
333bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
334bfcd1a73SHong Zhang   }
335bfcd1a73SHong Zhang #endif
336598bc09dSHong Zhang   PetscFunctionReturn(0);
337598bc09dSHong Zhang }
338598bc09dSHong Zhang 
3399123193aSHong Zhang #undef __FUNCT__
3409123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3419123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3429123193aSHong Zhang {
3439123193aSHong Zhang   PetscErrorCode ierr;
3449123193aSHong Zhang 
3459123193aSHong Zhang   PetscFunctionBegin;
3469123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
3479123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3489123193aSHong Zhang   }
3499123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3509123193aSHong Zhang   PetscFunctionReturn(0);
3519123193aSHong Zhang }
3529123193aSHong Zhang 
3534324174eSBarry Smith typedef struct {
3544324174eSBarry Smith   Mat         workB;
3554324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3564324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3574324174eSBarry Smith } MPIAIJ_MPIDense;
3584324174eSBarry Smith 
3597af9e4fdSHong Zhang #undef __FUNCT__
360*96b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
361*96b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3624324174eSBarry Smith {
3634324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3644324174eSBarry Smith   PetscErrorCode  ierr;
3654324174eSBarry Smith 
3664324174eSBarry Smith   PetscFunctionBegin;
3676bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3684324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3694324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3704324174eSBarry Smith   PetscFunctionReturn(0);
3714324174eSBarry Smith }
3724324174eSBarry Smith 
3739123193aSHong Zhang #undef __FUNCT__
3749123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
3759123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
3769123193aSHong Zhang {
3779123193aSHong Zhang   PetscErrorCode         ierr;
378f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
379d0f46423SBarry Smith   PetscInt               nz = aij->B->cmap->n;
380bf0cc555SLisandro Dalcin   PetscContainer         container;
3814324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
3824324174eSBarry Smith   VecScatter             ctx = aij->Mvctx;
3834324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
3844324174eSBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
385d0f46423SBarry Smith   PetscInt               m=A->rmap->n,n=B->cmap->n;
3869123193aSHong Zhang 
3879123193aSHong Zhang   PetscFunctionBegin;
388cb2480eaSBarry Smith   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
389d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
390cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
3911d567fd6SHong Zhang   ierr = MatMPIDenseSetPreallocation(*C,PETSC_NULL);CHKERRQ(ierr);
392cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
393cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3948cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense;
395f32d5d43SBarry Smith 
3964324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
397f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
398d0f46423SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
3994324174eSBarry Smith   /* Create work arrays needed */
400d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
401d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
4024324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
4034324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
4044324174eSBarry Smith 
405bf0cc555SLisandro Dalcin   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
406bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
407*96b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
408bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
409bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4109123193aSHong Zhang   PetscFunctionReturn(0);
4119123193aSHong Zhang }
4129123193aSHong Zhang 
4137af9e4fdSHong Zhang #undef __FUNCT__
4147af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
415f32d5d43SBarry Smith /*
4162636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4172636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
418f32d5d43SBarry Smith */
4194324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
420f32d5d43SBarry Smith {
421f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
422f32d5d43SBarry Smith   PetscErrorCode         ierr;
423f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
424f32d5d43SBarry Smith   VecScatter             ctx = aij->Mvctx;
425f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
426f32d5d43SBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
427f32d5d43SBarry Smith   PetscInt               i,j,k;
428aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
429aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
430f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
4317adad957SLisandro Dalcin   MPI_Comm               comm = ((PetscObject)A)->comm;
432d0f46423SBarry Smith   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
433f32d5d43SBarry Smith   MPI_Status             status;
4344324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
435bf0cc555SLisandro Dalcin   PetscContainer         container;
4364324174eSBarry Smith   Mat                    workB;
437f32d5d43SBarry Smith 
438f32d5d43SBarry Smith   PetscFunctionBegin;
439bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
44029825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
441bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4424324174eSBarry Smith 
4434324174eSBarry Smith   workB = *outworkB = contents->workB;
444cf1a0672SHong Zhang   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
445f32d5d43SBarry Smith   sindices  = to->indices;
446f32d5d43SBarry Smith   sstarts   = to->starts;
447f32d5d43SBarry Smith   sprocs    = to->procs;
4484324174eSBarry Smith   swaits    = contents->swaits;
4494324174eSBarry Smith   svalues   = contents->svalues;
450f32d5d43SBarry Smith 
451f32d5d43SBarry Smith   rindices  = from->indices;
452f32d5d43SBarry Smith   rstarts   = from->starts;
453f32d5d43SBarry Smith   rprocs    = from->procs;
4544324174eSBarry Smith   rwaits    = contents->rwaits;
4554324174eSBarry Smith   rvalues   = contents->rvalues;
456f32d5d43SBarry Smith 
457f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
458f32d5d43SBarry Smith   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
459f32d5d43SBarry Smith 
460f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
461f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
462f32d5d43SBarry Smith   }
463f32d5d43SBarry Smith 
464f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
465f32d5d43SBarry Smith     /* pack a message at a time */
466f32d5d43SBarry Smith     CHKMEMQ;
467f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
468f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
469f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
4702636ff26SBarry Smith       }
4712636ff26SBarry Smith     }
472f32d5d43SBarry Smith     CHKMEMQ;
473f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
474f32d5d43SBarry Smith   }
4752636ff26SBarry Smith 
476f32d5d43SBarry Smith   nrecvs = from->n;
477f32d5d43SBarry Smith   while (nrecvs) {
478f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
479f32d5d43SBarry Smith     nrecvs--;
480f32d5d43SBarry Smith     /* unpack a message at a time */
481f32d5d43SBarry Smith     CHKMEMQ;
482f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
483f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
484f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
4852636ff26SBarry Smith       }
4862636ff26SBarry Smith     }
487f32d5d43SBarry Smith     CHKMEMQ;
488f32d5d43SBarry Smith   }
489cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
490f32d5d43SBarry Smith 
491f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
492f32d5d43SBarry Smith   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
493f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
494f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
495f32d5d43SBarry Smith   PetscFunctionReturn(0);
496f32d5d43SBarry Smith }
497f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
498f32d5d43SBarry Smith 
4999123193aSHong Zhang #undef __FUNCT__
5009123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5019123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5029123193aSHong Zhang {
5039123193aSHong Zhang   PetscErrorCode       ierr;
504f32d5d43SBarry Smith   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
505f32d5d43SBarry Smith   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
506f32d5d43SBarry Smith   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
507f32d5d43SBarry Smith   Mat                  workB;
5089123193aSHong Zhang 
5099123193aSHong Zhang   PetscFunctionBegin;
5109123193aSHong Zhang 
511f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
512f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
513f32d5d43SBarry Smith 
514f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5154324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
516f32d5d43SBarry Smith 
517f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
518f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5199123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5209123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5219123193aSHong Zhang   PetscFunctionReturn(0);
5229123193aSHong Zhang }
523cf1a0672SHong Zhang 
5241634b032SHong Zhang #undef __FUNCT__
525cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
526cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
5271634b032SHong Zhang {
5281634b032SHong Zhang   PetscErrorCode     ierr;
529cf1a0672SHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
530cf1a0672SHong Zhang   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
531cf1a0672SHong Zhang   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
532cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
533cf1a0672SHong Zhang   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
534cf1a0672SHong Zhang   Mat_SeqAIJ         *p_loc,*p_oth;
535cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
53629825010SHong Zhang   PetscScalar        *pa_loc,*pa_oth,*pa,valtmp,*ca;
53729825010SHong Zhang   PetscInt           cm=C->rmap->n,anz,pnz;
538cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap=c->ptap;
53929825010SHong Zhang   PetscScalar        *apa_sparse=ptap->apa;
540cf1a0672SHong Zhang   PetscInt           *api,*apj,*apJ,i,j,k,row;
54129825010SHong Zhang   PetscInt           cstart=C->cmap->rstart;
54229825010SHong Zhang   PetscInt           cdnz,conz,k0,k1,nextp;
5431634b032SHong Zhang 
5441634b032SHong Zhang   PetscFunctionBegin;
545cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
546cf1a0672SHong Zhang   /*-----------------------------------------------------*/
547cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
548b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
549cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
550cf1a0672SHong Zhang 
551cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
552cf1a0672SHong Zhang   /*----------------------------------------------------------*/
553cf1a0672SHong Zhang   /* get data from symbolic products */
554cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
555cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
556cf1a0672SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
557cf1a0672SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
558cf1a0672SHong Zhang 
55929825010SHong Zhang   api = ptap->api;
56029825010SHong Zhang   apj = ptap->apj;
561cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
56229825010SHong Zhang     apJ = apj + api[i];
56329825010SHong Zhang 
564cf1a0672SHong Zhang     /* diagonal portion of A */
565cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
566cf1a0672SHong Zhang     adj = ad->j + adi[i];
567cf1a0672SHong Zhang     ada = ad->a + adi[i];
568cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
569cf1a0672SHong Zhang       row = adj[j];
570cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
571cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
572cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
57329825010SHong Zhang       /* perform sparse axpy */
574cf1a0672SHong Zhang       valtmp = ada[j];
57529825010SHong Zhang       nextp  = 0;
57629825010SHong Zhang       for (k=0; nextp<pnz; k++) {
57729825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
57829825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
57929825010SHong Zhang         }
5801634b032SHong Zhang       }
581cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
582cf1a0672SHong Zhang     }
5831634b032SHong Zhang 
584cf1a0672SHong Zhang     /* off-diagonal portion of A */
585cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
586cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
587cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
588cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
589cf1a0672SHong Zhang       row = aoj[j];
590cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
591cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
592cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
59329825010SHong Zhang       /* perform sparse axpy */
594cf1a0672SHong Zhang       valtmp = aoa[j];
59529825010SHong Zhang       nextp  = 0;
59629825010SHong Zhang       for (k=0; nextp<pnz; k++) {
59729825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
59829825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
59929825010SHong Zhang         }
600cf1a0672SHong Zhang       }
601cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
602cf1a0672SHong Zhang     }
6031634b032SHong Zhang 
604cf1a0672SHong Zhang     /* set values in C */
605cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
606cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6071634b032SHong Zhang 
608cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
609cf1a0672SHong Zhang     ca = coa + co->i[i];
610cf1a0672SHong Zhang     k  = 0;
611cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++){
612cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
61329825010SHong Zhang       ca[k0]      = apa_sparse[k];
61429825010SHong Zhang       apa_sparse[k] = 0.0;
615cf1a0672SHong Zhang       k++;
616cf1a0672SHong Zhang     }
6171634b032SHong Zhang 
618cf1a0672SHong Zhang     /* diagonal part of C */
619cf1a0672SHong Zhang     ca = cda + cd->i[i];
620cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++){
62129825010SHong Zhang       ca[k1]      = apa_sparse[k];
62229825010SHong Zhang       apa_sparse[k] = 0.0;
623cf1a0672SHong Zhang       k++;
624cf1a0672SHong Zhang     }
625cf1a0672SHong Zhang 
626cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
627cf1a0672SHong Zhang     ca = coa + co->i[i];
628cf1a0672SHong Zhang     for (; k0<conz; k0++){
62929825010SHong Zhang       ca[k0]      = apa_sparse[k];
63029825010SHong Zhang       apa_sparse[k] = 0.0;
631cf1a0672SHong Zhang       k++;
632cf1a0672SHong Zhang     }
633cf1a0672SHong Zhang   }
634cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
635cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
636cf1a0672SHong Zhang   PetscFunctionReturn(0);
637cf1a0672SHong Zhang }
638cf1a0672SHong Zhang 
639cf1a0672SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */
640cf1a0672SHong Zhang #undef __FUNCT__
641cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
642cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C)
643cf1a0672SHong Zhang {
644cf1a0672SHong Zhang   PetscErrorCode       ierr;
645cf1a0672SHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
646cf1a0672SHong Zhang   Mat                  Cmpi;
647cf1a0672SHong Zhang   Mat_PtAPMPI          *ptap;
648cf1a0672SHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
649cf1a0672SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
650cf1a0672SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
651cf1a0672SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
652cf1a0672SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
653f84c536bSHong Zhang   PetscInt             i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
654cf1a0672SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
655cf1a0672SHong Zhang   PetscInt             nlnk_max,armax,prmax;
656cf1a0672SHong Zhang   PetscReal            afill;
657cf1a0672SHong Zhang   PetscScalar          *apa;
658cf1a0672SHong Zhang 
659cf1a0672SHong Zhang   PetscFunctionBegin;
660cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
661cf1a0672SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
662cf1a0672SHong Zhang 
663cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
664b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
66519f0e57aSHong Zhang 
666cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
667cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
668cf1a0672SHong Zhang 
669cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
670cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
671cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
672cf1a0672SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
673cf1a0672SHong Zhang 
674cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
675cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
676cf1a0672SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
677cf1a0672SHong Zhang   ptap->api = api;
678cf1a0672SHong Zhang   api[0]    = 0;
679cf1a0672SHong Zhang 
680cf1a0672SHong Zhang   /* create and initialize a linked list */
681cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
682cf1a0672SHong Zhang   prmax = PetscMax(p_loc->rmax,p_oth->rmax);
683cf1a0672SHong Zhang   nlnk_max = armax*prmax;
684cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
685f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
686cf1a0672SHong Zhang 
687cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
688cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
689cf1a0672SHong Zhang   current_space = free_space;
690cf1a0672SHong Zhang 
691cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
692cf1a0672SHong Zhang   for (i=0; i<am; i++) {
693cf1a0672SHong Zhang     apnz = 0;
694cf1a0672SHong Zhang     /* diagonal portion of A */
695cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
696cf1a0672SHong Zhang     for (j=0; j<nzi; j++){
697cf1a0672SHong Zhang       row = *adj++;
698cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
699cf1a0672SHong Zhang       Jptr  = pj_loc + pi_loc[row];
700cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
701f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
702cf1a0672SHong Zhang     }
703cf1a0672SHong Zhang     /* off-diagonal portion of A */
704cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
705cf1a0672SHong Zhang     for (j=0; j<nzi; j++){
706cf1a0672SHong Zhang       row = *aoj++;
707cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
708cf1a0672SHong Zhang       Jptr  = pj_oth + pi_oth[row];
709f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
710cf1a0672SHong Zhang     }
711cf1a0672SHong Zhang 
712f84c536bSHong Zhang     apnz     = *lnk;
713cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
714cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
715cf1a0672SHong Zhang 
716cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
717cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
718cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
719cf1a0672SHong Zhang       nspacedouble++;
720cf1a0672SHong Zhang     }
721cf1a0672SHong Zhang 
722cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
723f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
724cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
725cf1a0672SHong Zhang     current_space->array           += apnz;
726cf1a0672SHong Zhang     current_space->local_used      += apnz;
727cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
728cf1a0672SHong Zhang   }
729cf1a0672SHong Zhang 
730cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
731cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
732cf1a0672SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
733cf1a0672SHong Zhang   apj  = ptap->apj;
734cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
735f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
736cf1a0672SHong Zhang 
737cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
738cf1a0672SHong Zhang   /*----------------------------------------------------*/
739cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
740cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
741cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
742cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
743cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
744cf1a0672SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
745cf1a0672SHong Zhang 
74629825010SHong Zhang   /* malloc apa for assembly Cmpi */
747cf1a0672SHong Zhang   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
748cf1a0672SHong Zhang   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
74929825010SHong Zhang   ptap->apa = apa;
750cf1a0672SHong Zhang   for (i=0; i<am; i++){
751cf1a0672SHong Zhang     row  = i + rstart;
752cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
753cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
754cf1a0672SHong Zhang     apj += apnz;
755cf1a0672SHong Zhang   }
756cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
757cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
758cf1a0672SHong Zhang 
759cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
760cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
761cf1a0672SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
762cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
763cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
764cf1a0672SHong Zhang 
765cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
766cf1a0672SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
767cf1a0672SHong Zhang   c->ptap  = ptap;
768cf1a0672SHong Zhang 
769cf1a0672SHong Zhang   *C = Cmpi;
770cf1a0672SHong Zhang 
771cf1a0672SHong Zhang   /* set MatInfo */
772cf1a0672SHong Zhang   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
773cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
774cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
775cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
776cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
777cf1a0672SHong Zhang 
778cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
779cf1a0672SHong Zhang   if (api[am]) {
780cf1a0672SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
781cf1a0672SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
782cf1a0672SHong Zhang   } else {
783cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
784cf1a0672SHong Zhang   }
785cf1a0672SHong Zhang #endif
7861634b032SHong Zhang   PetscFunctionReturn(0);
7871634b032SHong Zhang }
788f96d37fcSHong Zhang 
789f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
790f96d37fcSHong Zhang #undef __FUNCT__
791f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
792c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
793f96d37fcSHong Zhang {
794f96d37fcSHong Zhang   PetscErrorCode ierr;
795d6ab1dc0SHong Zhang   PetscBool      scalable=PETSC_FALSE;
796f96d37fcSHong Zhang 
797f96d37fcSHong Zhang   PetscFunctionBegin;
798f96d37fcSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
799d6ab1dc0SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
800d6ab1dc0SHong Zhang       ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
801d6ab1dc0SHong Zhang       if  (scalable){
802d6ab1dc0SHong Zhang         ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr);
803d6ab1dc0SHong Zhang       } else {
804c5af039cSHong Zhang         ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
805f96d37fcSHong Zhang       }
806d6ab1dc0SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
807d6ab1dc0SHong Zhang   }
808d6ab1dc0SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
809f96d37fcSHong Zhang   PetscFunctionReturn(0);
810f96d37fcSHong Zhang }
811f96d37fcSHong Zhang 
812f96d37fcSHong Zhang #undef __FUNCT__
813f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
814c5af039cSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
815f96d37fcSHong Zhang {
816c5af039cSHong Zhang   PetscErrorCode       ierr;
817c5af039cSHong Zhang   Mat_Merge_SeqsToMPI  *merge;
818497f5370SHong Zhang   Mat_MPIAIJ           *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
819c5af039cSHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
820c5af039cSHong Zhang   Mat_PtAPMPI          *ptap;
821d6ab1dc0SHong Zhang   PetscInt             *adj,*aJ;
822497f5370SHong Zhang   PetscInt             i,j,k,anz,pnz,row,*cj;
823e745005bSHong Zhang   MatScalar            *ada,*aval,*ca,valtmp;
824c5af039cSHong Zhang   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
825c5af039cSHong Zhang   MPI_Comm             comm=((PetscObject)C)->comm;
826c5af039cSHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
827c5af039cSHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
828c5af039cSHong Zhang   PetscInt             **buf_ri,**buf_rj;
829c5af039cSHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
830c5af039cSHong Zhang   MPI_Request          *s_waits,*r_waits;
831c5af039cSHong Zhang   MPI_Status           *status;
832c5af039cSHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
833d6ab1dc0SHong Zhang   PetscInt             *ai,*aj,*coi,*coj;
834497f5370SHong Zhang   PetscInt             *poJ=po->j,*pdJ=pd->j;
835c5af039cSHong Zhang   Mat                  A_loc;
836c5af039cSHong Zhang   Mat_SeqAIJ           *a_loc;
837f96d37fcSHong Zhang 
838f96d37fcSHong Zhang   PetscFunctionBegin;
839c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
840c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
841c5af039cSHong Zhang 
842c5af039cSHong Zhang   ptap  = c->ptap;
843c5af039cSHong Zhang   merge = ptap->merge;
844c5af039cSHong Zhang 
845c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
846c5af039cSHong Zhang   /*--------------------------------------------------------------*/
847c5af039cSHong Zhang   /* get data from symbolic products */
848c5af039cSHong Zhang   coi = merge->coi; coj = merge->coj;
849c5af039cSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
850c5af039cSHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
851c5af039cSHong Zhang 
852c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
853c5af039cSHong Zhang   owners = merge->rowmap->range;
854c5af039cSHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
855c5af039cSHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
856c5af039cSHong Zhang 
857c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
858c5af039cSHong Zhang   A_loc = ptap->A_loc;
859c5af039cSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
860c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
861d6ab1dc0SHong Zhang   ai   = a_loc->i;
862d6ab1dc0SHong Zhang   aj   = a_loc->j;
863c5af039cSHong Zhang 
864e745005bSHong Zhang   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */
865e745005bSHong Zhang   ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
866c5af039cSHong Zhang 
867c5af039cSHong Zhang     for (i=0; i<am; i++) {
868e745005bSHong Zhang       /* 2-a) put A[i,:] to dense array aval */
869d6ab1dc0SHong Zhang       anz = ai[i+1] - ai[i];
870d6ab1dc0SHong Zhang       adj = aj + ai[i];
871d6ab1dc0SHong Zhang       ada = a_loc->a + ai[i];
872c5af039cSHong Zhang       for (j=0; j<anz; j++){
873e745005bSHong Zhang         aval[adj[j]] = ada[j];
874c5af039cSHong Zhang       }
875c5af039cSHong Zhang 
876c5af039cSHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
877c5af039cSHong Zhang       /*--------------------------------------------------------------*/
878c5af039cSHong Zhang       /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
879c5af039cSHong Zhang       pnz = po->i[i+1] - po->i[i];
880c5af039cSHong Zhang       poJ = po->j + po->i[i];
881c5af039cSHong Zhang       pA  = po->a + po->i[i];
882c5af039cSHong Zhang       for (j=0; j<pnz; j++){
883c5af039cSHong Zhang         row = poJ[j];
884c5af039cSHong Zhang         cnz = coi[row+1] - coi[row];
885c5af039cSHong Zhang         cj  = coj + coi[row];
886c5af039cSHong Zhang         ca  = coa + coi[row];
887c5af039cSHong Zhang         /* perform dense axpy */
888c5af039cSHong Zhang         valtmp = pA[j];
889c5af039cSHong Zhang         for (k=0; k<cnz; k++) {
890e745005bSHong Zhang           ca[k] += valtmp*aval[cj[k]];
891c5af039cSHong Zhang         }
892c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
893c5af039cSHong Zhang       }
894c5af039cSHong Zhang 
895c5af039cSHong Zhang       /* put the value into Cd (diagonal part) */
896c5af039cSHong Zhang       pnz = pd->i[i+1] - pd->i[i];
897c5af039cSHong Zhang       pdJ = pd->j + pd->i[i];
898c5af039cSHong Zhang       pA  = pd->a + pd->i[i];
899c5af039cSHong Zhang       for (j=0; j<pnz; j++){
900c5af039cSHong Zhang         row = pdJ[j];
901c5af039cSHong Zhang         cnz = bi[row+1] - bi[row];
902c5af039cSHong Zhang         cj  = bj + bi[row];
903c5af039cSHong Zhang         ca  = ba + bi[row];
904c5af039cSHong Zhang         /* perform dense axpy */
905c5af039cSHong Zhang         valtmp = pA[j];
906c5af039cSHong Zhang         for (k=0; k<cnz; k++) {
907e745005bSHong Zhang           ca[k] += valtmp*aval[cj[k]];
908c5af039cSHong Zhang         }
909c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
910c5af039cSHong Zhang       }
911c5af039cSHong Zhang 
912d6ab1dc0SHong Zhang       /* zero the current row of Pt*A */
913d6ab1dc0SHong Zhang       aJ = aj + ai[i];
914e745005bSHong Zhang       for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
915c5af039cSHong Zhang     }
916c5af039cSHong Zhang 
917c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
918c5af039cSHong Zhang   /*------------------------------------*/
919c5af039cSHong Zhang   buf_ri = merge->buf_ri;
920c5af039cSHong Zhang   buf_rj = merge->buf_rj;
921c5af039cSHong Zhang   len_s  = merge->len_s;
922c5af039cSHong Zhang   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
923c5af039cSHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
924c5af039cSHong Zhang 
925c5af039cSHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
926c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++){
927c5af039cSHong Zhang     if (!len_s[proc]) continue;
928c5af039cSHong Zhang     i = merge->owners_co[proc];
929c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
930c5af039cSHong Zhang     k++;
931c5af039cSHong Zhang   }
932c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
933c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
934c5af039cSHong Zhang 
935c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
936c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
937c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
938c5af039cSHong Zhang 
939c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
940c5af039cSHong Zhang   /*----------------------------------------------------*/
941c5af039cSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
942c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++){
943c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
944c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
945c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
946c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
947c5af039cSHong Zhang   }
948c5af039cSHong Zhang 
949c5af039cSHong Zhang   for (i=0; i<cm; i++) {
950c5af039cSHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
951c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
952c5af039cSHong Zhang     ba_i = ba + bi[i];
953c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
954c5af039cSHong Zhang     /* add received vals into ba */
955c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
956c5af039cSHong Zhang       /* i-th row */
957c5af039cSHong Zhang       if (i == *nextrow[k]) {
958c5af039cSHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
959c5af039cSHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
960c5af039cSHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
961c5af039cSHong Zhang         nextcj = 0;
962c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++){
963c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
964c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
965c5af039cSHong Zhang           }
966c5af039cSHong Zhang         }
967c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
968c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
969c5af039cSHong Zhang       }
970c5af039cSHong Zhang     }
971c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
972c5af039cSHong Zhang   }
973c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
974c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
975c5af039cSHong Zhang 
976c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
977c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
978c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
979c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
980e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
981f96d37fcSHong Zhang   PetscFunctionReturn(0);
982f96d37fcSHong Zhang }
983f96d37fcSHong Zhang 
984c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
985f96d37fcSHong Zhang #undef __FUNCT__
986f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
987f96d37fcSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
988f96d37fcSHong Zhang {
989f96d37fcSHong Zhang   PetscErrorCode       ierr;
9904e938277SHong Zhang   Mat                  Cmpi,A_loc,POt,PDt;
991f96d37fcSHong Zhang   Mat_PtAPMPI          *ptap;
992f96d37fcSHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
993f96d37fcSHong Zhang   Mat_MPIAIJ           *p=(Mat_MPIAIJ*)P->data,*c;
994f96d37fcSHong Zhang   PetscInt             *pdti,*pdtj,*poti,*potj,*ptJ;
995f96d37fcSHong Zhang   PetscInt             nnz;
9964e938277SHong Zhang   PetscInt             *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
997497f5370SHong Zhang   PetscInt             am=A->rmap->n,pn=P->cmap->n;
998f96d37fcSHong Zhang   PetscBT              lnkbt;
999f96d37fcSHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
1000f96d37fcSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1001f96d37fcSHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
1002f96d37fcSHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
1003f96d37fcSHong Zhang   PetscInt             nzi,*bi,*bj;
1004f96d37fcSHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1005f96d37fcSHong Zhang   MPI_Request          *swaits,*rwaits;
1006f96d37fcSHong Zhang   MPI_Status           *sstatus,rstatus;
1007f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI  *merge;
1008d6ab1dc0SHong Zhang   PetscInt             *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1009f96d37fcSHong Zhang   PetscReal            afill=1.0,afill_tmp;
1010c36aecf5SHong Zhang   PetscInt             rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1011f96d37fcSHong Zhang   PetscScalar          *vals;
10124e938277SHong Zhang   Mat_SeqAIJ           *a_loc, *pdt,*pot;
1013f96d37fcSHong Zhang 
1014f96d37fcSHong Zhang   PetscFunctionBegin;
1015c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1016c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
1017c5af039cSHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1018c5af039cSHong Zhang   }
1019c5af039cSHong Zhang 
1020f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1021f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1022f96d37fcSHong Zhang 
1023f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1024f96d37fcSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1025f96d37fcSHong Zhang 
1026f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1027f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1028c5af039cSHong Zhang   ptap->A_loc = A_loc;
10291c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1030d6ab1dc0SHong Zhang   ai   = a_loc->i;
1031d6ab1dc0SHong Zhang   aj   = a_loc->j;
1032f96d37fcSHong Zhang 
1033f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1034f96d37fcSHong Zhang   /*----------------------------------------------------*/
10354e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
10364e938277SHong Zhang   pdt = (Mat_SeqAIJ*)PDt->data;
10374e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
10384e938277SHong Zhang 
10394e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
10404e938277SHong Zhang   pot = (Mat_SeqAIJ*)POt->data;
10414e938277SHong Zhang   poti = pot->i; potj = pot->j;
1042f96d37fcSHong Zhang 
1043f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1044f96d37fcSHong Zhang   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1045f96d37fcSHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1046f96d37fcSHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1047f96d37fcSHong Zhang   coi[0] = 0;
1048f96d37fcSHong Zhang 
1049f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1050d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1051f96d37fcSHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1052f96d37fcSHong Zhang   current_space = free_space;
105319f0e57aSHong Zhang 
1054f96d37fcSHong Zhang   /* create and initialize a linked list */
10554e938277SHong Zhang   i = PetscMax(pdt->rmax,pot->rmax);
1056c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
10574e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
10584e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1059f96d37fcSHong Zhang 
1060f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1061f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1062f96d37fcSHong Zhang     ptJ = potj + poti[i];
1063f96d37fcSHong Zhang     for (j=0; j<pnz; j++){
1064f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1065d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1066d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1067f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1068d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1069f96d37fcSHong Zhang     }
10704e938277SHong Zhang     nnz = lnk[0];
1071f96d37fcSHong Zhang 
1072f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1073f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1074f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1075f96d37fcSHong Zhang       nspacedouble++;
1076f96d37fcSHong Zhang     }
1077f96d37fcSHong Zhang 
1078f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
10794e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1080f96d37fcSHong Zhang     current_space->array           += nnz;
1081f96d37fcSHong Zhang     current_space->local_used      += nnz;
1082f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
1083f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1084f96d37fcSHong Zhang   }
1085f96d37fcSHong Zhang 
1086f96d37fcSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1087f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1088d6ab1dc0SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]);
1089f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1090f96d37fcSHong Zhang 
1091f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1092f96d37fcSHong Zhang   /*----------------------------------------------*/
1093f96d37fcSHong Zhang   /* determine row ownership */
1094f96d37fcSHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1095f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1096f96d37fcSHong Zhang   merge->rowmap->n = pn;
1097f96d37fcSHong Zhang   merge->rowmap->bs = 1;
1098f96d37fcSHong Zhang   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1099f96d37fcSHong Zhang   owners = merge->rowmap->range;
1100f96d37fcSHong Zhang 
1101f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
1102f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1103f96d37fcSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1104f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1105f96d37fcSHong Zhang   len_s = merge->len_s;
1106f96d37fcSHong Zhang   merge->nsend = 0;
1107f96d37fcSHong Zhang 
1108f96d37fcSHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1109f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1110f96d37fcSHong Zhang 
1111f96d37fcSHong Zhang   proc = 0;
1112f96d37fcSHong Zhang   for (i=0; i<pon; i++){
1113f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1114f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1115f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1116f96d37fcSHong Zhang   }
1117f96d37fcSHong Zhang 
1118f96d37fcSHong Zhang   len   = 0;  /* max length of buf_si[] */
1119f96d37fcSHong Zhang   owners_co[0] = 0;
1120f96d37fcSHong Zhang   for (proc=0; proc<size; proc++){
1121f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1122f96d37fcSHong Zhang     if (len_si[proc]){
1123f96d37fcSHong Zhang       merge->nsend++;
1124f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1125f96d37fcSHong Zhang       len += len_si[proc];
1126f96d37fcSHong Zhang     }
1127f96d37fcSHong Zhang   }
1128f96d37fcSHong Zhang 
1129f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
1130f96d37fcSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1131f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1132f96d37fcSHong Zhang 
1133f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1134f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1135f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1136f96d37fcSHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1137f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++){
1138f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1139f96d37fcSHong Zhang     i = owners_co[proc];
1140f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1141f96d37fcSHong Zhang     k++;
1142f96d37fcSHong Zhang   }
1143f96d37fcSHong Zhang 
1144f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1145f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1146f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++){
1147c280ed6aSJed Brown     PetscMPIInt icompleted;
1148c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1149f96d37fcSHong Zhang   }
1150f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1151f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1152f96d37fcSHong Zhang 
1153f96d37fcSHong Zhang   /* send and recv coi */
1154f96d37fcSHong Zhang   /*-------------------*/
1155f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1156f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1157f96d37fcSHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1158f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1159f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++){
1160f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1161f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1162f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1163f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1164f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1165f96d37fcSHong Zhang     */
1166f96d37fcSHong Zhang     /*-------------------------------------------*/
1167f96d37fcSHong Zhang     nrows = len_si[proc]/2 - 1;
1168f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1169f96d37fcSHong Zhang     buf_si[0]   = nrows;
1170f96d37fcSHong Zhang     buf_si_i[0] = 0;
1171f96d37fcSHong Zhang     nrows = 0;
1172f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
1173f96d37fcSHong Zhang       nzi = coi[i+1] - coi[i];
1174f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1175f96d37fcSHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
1176f96d37fcSHong Zhang       nrows++;
1177f96d37fcSHong Zhang     }
1178f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1179f96d37fcSHong Zhang     k++;
1180f96d37fcSHong Zhang     buf_si += len_si[proc];
1181f96d37fcSHong Zhang   }
1182f96d37fcSHong Zhang   i = merge->nrecv;
1183f96d37fcSHong Zhang   while (i--) {
1184c280ed6aSJed Brown     PetscMPIInt icompleted;
1185c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1186f96d37fcSHong Zhang   }
1187f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1188f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1189f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1190f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1191f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1192f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1193f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1194f96d37fcSHong Zhang 
1195f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1196f96d37fcSHong Zhang   /*------------------------------------------*/
1197f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1198f96d37fcSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1199f96d37fcSHong Zhang   bi[0] = 0;
1200f96d37fcSHong Zhang 
1201c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1202d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1203f96d37fcSHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1204f96d37fcSHong Zhang   current_space = free_space;
1205f96d37fcSHong Zhang 
1206f96d37fcSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1207f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++){
1208f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1209f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1210f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1211f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
1212f96d37fcSHong Zhang   }
1213f96d37fcSHong Zhang 
12141c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1215f96d37fcSHong Zhang   rmax = 0;
1216f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1217f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1218f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1219f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1220f96d37fcSHong Zhang     for (j=0; j<pnz; j++){
1221f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1222d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1223d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1224f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1225d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1226f96d37fcSHong Zhang     }
12274e938277SHong Zhang 
1228f96d37fcSHong Zhang     /* add received col data into lnk */
1229f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
1230f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1231f96d37fcSHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
1232f96d37fcSHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
12334e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1234f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1235f96d37fcSHong Zhang       }
1236f96d37fcSHong Zhang     }
12374e938277SHong Zhang     nnz = lnk[0];
1238f96d37fcSHong Zhang 
1239f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1240f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1241f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1242f96d37fcSHong Zhang       nspacedouble++;
1243f96d37fcSHong Zhang     }
1244f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
12454e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1246f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1247f96d37fcSHong Zhang     current_space->array           += nnz;
1248f96d37fcSHong Zhang     current_space->local_used      += nnz;
1249f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
1250f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1251f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1252f96d37fcSHong Zhang   }
1253f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1254f96d37fcSHong Zhang 
1255f96d37fcSHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1256f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1257d6ab1dc0SHong Zhang   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]);
1258f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1259d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
12604e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
12614e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1262f96d37fcSHong Zhang 
12631c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
12641c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
12651c7d5954SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
12661c7d5954SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1267f96d37fcSHong Zhang 
1268f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1269f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1270f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1271f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1272f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1273f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1274f96d37fcSHong Zhang   for (i=0; i<pn; i++){
1275f96d37fcSHong Zhang     row = i + rstart;
1276f96d37fcSHong Zhang     nnz = bi[i+1] - bi[i];
1277f96d37fcSHong Zhang     Jptr = bj + bi[i];
1278f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1279f96d37fcSHong Zhang   }
1280f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1281f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1282f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1283f96d37fcSHong Zhang 
1284f96d37fcSHong Zhang   merge->bi            = bi;
1285f96d37fcSHong Zhang   merge->bj            = bj;
1286f96d37fcSHong Zhang   merge->coi           = coi;
1287f96d37fcSHong Zhang   merge->coj           = coj;
1288f96d37fcSHong Zhang   merge->buf_ri        = buf_ri;
1289f96d37fcSHong Zhang   merge->buf_rj        = buf_rj;
1290f96d37fcSHong Zhang   merge->owners_co     = owners_co;
1291f96d37fcSHong Zhang   merge->destroy       = Cmpi->ops->destroy;
1292f96d37fcSHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
1293f96d37fcSHong Zhang 
1294d6ab1dc0SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1295c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1296f96d37fcSHong Zhang 
1297f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1298f96d37fcSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
1299f96d37fcSHong Zhang   c->ptap        = ptap;
1300c5af039cSHong Zhang   ptap->api      = PETSC_NULL;
1301c5af039cSHong Zhang   ptap->apj      = PETSC_NULL;
1302f96d37fcSHong Zhang   ptap->merge    = merge;
1303d6ab1dc0SHong Zhang   ptap->rmax     = rmax;
1304d6ab1dc0SHong Zhang 
1305d6ab1dc0SHong Zhang   *C = Cmpi;
1306d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1307d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
1308d6ab1dc0SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1309d6ab1dc0SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1310d6ab1dc0SHong Zhang   } else {
1311d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1312d6ab1dc0SHong Zhang   }
1313d6ab1dc0SHong Zhang #endif
1314d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1315d6ab1dc0SHong Zhang }
1316d6ab1dc0SHong Zhang 
1317d6ab1dc0SHong Zhang #undef __FUNCT__
1318d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
1319d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C)
1320d6ab1dc0SHong Zhang {
1321d6ab1dc0SHong Zhang   PetscErrorCode       ierr;
1322d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
1323d6ab1dc0SHong Zhang   Mat_MPIAIJ           *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1324d6ab1dc0SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1325d6ab1dc0SHong Zhang   Mat_PtAPMPI          *ptap;
1326e745005bSHong Zhang   PetscInt             *adj;
1327e745005bSHong Zhang   PetscInt             i,j,k,anz,pnz,row,*cj,nexta;
1328e745005bSHong Zhang   MatScalar            *ada,*ca,valtmp;
1329d6ab1dc0SHong Zhang   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1330d6ab1dc0SHong Zhang   MPI_Comm             comm=((PetscObject)C)->comm;
1331d6ab1dc0SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
1332d6ab1dc0SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1333d6ab1dc0SHong Zhang   PetscInt             **buf_ri,**buf_rj;
1334d6ab1dc0SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
1335d6ab1dc0SHong Zhang   MPI_Request          *s_waits,*r_waits;
1336d6ab1dc0SHong Zhang   MPI_Status           *status;
1337d6ab1dc0SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
1338d6ab1dc0SHong Zhang   PetscInt             *ai,*aj,*coi,*coj;
1339d6ab1dc0SHong Zhang   PetscInt             *poJ=po->j,*pdJ=pd->j;
1340d6ab1dc0SHong Zhang   Mat                  A_loc;
1341d6ab1dc0SHong Zhang   Mat_SeqAIJ           *a_loc;
1342d6ab1dc0SHong Zhang 
1343d6ab1dc0SHong Zhang   PetscFunctionBegin;
1344d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1345d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1346d6ab1dc0SHong Zhang 
1347d6ab1dc0SHong Zhang   ptap  = c->ptap;
1348d6ab1dc0SHong Zhang   merge = ptap->merge;
1349d6ab1dc0SHong Zhang 
1350e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1351e745005bSHong Zhang   /*------------------------------------------*/
1352d6ab1dc0SHong Zhang   /* get data from symbolic products */
1353d6ab1dc0SHong Zhang   coi = merge->coi; coj = merge->coj;
1354d6ab1dc0SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1355d6ab1dc0SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1356d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1357d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1358d6ab1dc0SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1359d6ab1dc0SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1360d6ab1dc0SHong Zhang 
1361d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1362d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1363d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1364d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1365d6ab1dc0SHong Zhang   ai   = a_loc->i;
1366d6ab1dc0SHong Zhang   aj   = a_loc->j;
1367d6ab1dc0SHong Zhang 
1368d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1369e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
1370d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1371d6ab1dc0SHong Zhang     adj = aj + ai[i];
1372d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1373d6ab1dc0SHong Zhang 
1374d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1375e745005bSHong Zhang     /*-------------------------------------------------------------*/
1376d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1377d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1378d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1379d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1380d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++){
1381d6ab1dc0SHong Zhang       row = poJ[j];
1382d6ab1dc0SHong Zhang       cnz = coi[row+1] - coi[row];
1383d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1384d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1385e745005bSHong Zhang       /* perform sparse axpy */
1386e745005bSHong Zhang       nexta  = 0;
1387d6ab1dc0SHong Zhang       valtmp = pA[j];
1388e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1389e745005bSHong Zhang         if (cj[k] == adj[nexta]){
1390e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1391e745005bSHong Zhang           nexta++;
1392d6ab1dc0SHong Zhang         }
1393e745005bSHong Zhang       }
1394e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1395d6ab1dc0SHong Zhang     }
1396d6ab1dc0SHong Zhang 
1397d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1398d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1399d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1400d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1401d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++){
1402d6ab1dc0SHong Zhang       row = pdJ[j];
1403d6ab1dc0SHong Zhang       cnz = bi[row+1] - bi[row];
1404d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1405d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1406e745005bSHong Zhang       /* perform sparse axpy */
1407e745005bSHong Zhang       nexta  = 0;
1408d6ab1dc0SHong Zhang       valtmp = pA[j];
1409e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1410e745005bSHong Zhang         if (cj[k] == adj[nexta]){
1411e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1412e745005bSHong Zhang           nexta++;
1413d6ab1dc0SHong Zhang         }
1414e745005bSHong Zhang       }
1415e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1416d6ab1dc0SHong Zhang     }
1417d6ab1dc0SHong Zhang 
1418d6ab1dc0SHong Zhang   }
1419d6ab1dc0SHong Zhang 
1420d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1421d6ab1dc0SHong Zhang   /*------------------------------------*/
1422d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1423d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1424d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1425d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1426d6ab1dc0SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1427d6ab1dc0SHong Zhang 
1428d6ab1dc0SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1429d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++){
1430d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1431d6ab1dc0SHong Zhang     i = merge->owners_co[proc];
1432d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1433d6ab1dc0SHong Zhang     k++;
1434d6ab1dc0SHong Zhang   }
1435d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1436d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1437d6ab1dc0SHong Zhang 
1438d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1439d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1440d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1441d6ab1dc0SHong Zhang 
1442d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1443d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1444d6ab1dc0SHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1445d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++){
1446e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1447d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1448d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1449d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
1450d6ab1dc0SHong Zhang   }
1451d6ab1dc0SHong Zhang 
1452d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1453d6ab1dc0SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
1454d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1455d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1456d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1457d6ab1dc0SHong Zhang     /* add received vals into ba */
1458d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
1459d6ab1dc0SHong Zhang       /* i-th row */
1460d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1461d6ab1dc0SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
1462d6ab1dc0SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
1463d6ab1dc0SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
1464d6ab1dc0SHong Zhang         nextcj = 0;
1465d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++){
1466d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
1467d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1468d6ab1dc0SHong Zhang           }
1469d6ab1dc0SHong Zhang         }
1470d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1471d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1472d6ab1dc0SHong Zhang       }
1473d6ab1dc0SHong Zhang     }
1474d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1475d6ab1dc0SHong Zhang   }
1476d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1477d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1478d6ab1dc0SHong Zhang 
1479d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1480d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1481d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1482d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1483d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1484d6ab1dc0SHong Zhang }
1485d6ab1dc0SHong Zhang 
1486d6ab1dc0SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1487d6ab1dc0SHong Zhang #undef __FUNCT__
1488d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
1489d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C)
1490d6ab1dc0SHong Zhang {
1491d6ab1dc0SHong Zhang   PetscErrorCode       ierr;
1492d6ab1dc0SHong Zhang   Mat                  Cmpi,A_loc,POt,PDt;
1493d6ab1dc0SHong Zhang   Mat_PtAPMPI          *ptap;
1494d6ab1dc0SHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
1495d6ab1dc0SHong Zhang   Mat_MPIAIJ           *p=(Mat_MPIAIJ*)P->data,*c;
1496d6ab1dc0SHong Zhang   PetscInt             *pdti,*pdtj,*poti,*potj,*ptJ;
1497d6ab1dc0SHong Zhang   PetscInt             nnz;
1498d6ab1dc0SHong Zhang   PetscInt             *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1499d6ab1dc0SHong Zhang   PetscInt             am=A->rmap->n,pn=P->cmap->n;
1500d6ab1dc0SHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
1501d6ab1dc0SHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1502d6ab1dc0SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
1503d6ab1dc0SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
1504d6ab1dc0SHong Zhang   PetscInt             nzi,*bi,*bj;
1505d6ab1dc0SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1506d6ab1dc0SHong Zhang   MPI_Request          *swaits,*rwaits;
1507d6ab1dc0SHong Zhang   MPI_Status           *sstatus,rstatus;
1508d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
1509d6ab1dc0SHong Zhang   PetscInt             *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1510d6ab1dc0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
1511c36aecf5SHong Zhang   PetscInt             rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1512d6ab1dc0SHong Zhang   PetscScalar          *vals;
1513d6ab1dc0SHong Zhang   Mat_SeqAIJ           *a_loc, *pdt,*pot;
1514d6ab1dc0SHong Zhang 
1515d6ab1dc0SHong Zhang   PetscFunctionBegin;
1516d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1517d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
1518d6ab1dc0SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1519d6ab1dc0SHong Zhang   }
1520d6ab1dc0SHong Zhang 
1521d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1522d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1523d6ab1dc0SHong Zhang 
1524d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1525d6ab1dc0SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1526d6ab1dc0SHong Zhang 
1527d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1528d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1530d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1531d6ab1dc0SHong Zhang   ai   = a_loc->i;
1532d6ab1dc0SHong Zhang   aj   = a_loc->j;
1533d6ab1dc0SHong Zhang 
1534d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1535d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1536d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1537d6ab1dc0SHong Zhang   pdt = (Mat_SeqAIJ*)PDt->data;
1538d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1539d6ab1dc0SHong Zhang 
1540d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1541d6ab1dc0SHong Zhang   pot = (Mat_SeqAIJ*)POt->data;
1542d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1543d6ab1dc0SHong Zhang 
1544d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1545d6ab1dc0SHong Zhang   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1546d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1547d6ab1dc0SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1548d6ab1dc0SHong Zhang   coi[0] = 0;
1549d6ab1dc0SHong Zhang 
1550d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1551d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1552d6ab1dc0SHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1553d6ab1dc0SHong Zhang   current_space = free_space;
155419f0e57aSHong Zhang 
1555d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1556d6ab1dc0SHong Zhang   i = PetscMax(pdt->rmax,pot->rmax);
1557c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1558d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1559d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1560d6ab1dc0SHong Zhang 
1561d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1562d6ab1dc0SHong Zhang     nnz = 0;
1563d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1564d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1565d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++){
1566d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1567d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1568d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1569d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1570d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1571d6ab1dc0SHong Zhang     }
1572d6ab1dc0SHong Zhang     nnz = lnk[0];
1573d6ab1dc0SHong Zhang 
1574d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1575d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1576d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1577d6ab1dc0SHong Zhang       nspacedouble++;
1578d6ab1dc0SHong Zhang     }
1579d6ab1dc0SHong Zhang 
1580d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1581d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1582d6ab1dc0SHong Zhang     current_space->array           += nnz;
1583d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1584d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
1585d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1586d6ab1dc0SHong Zhang   }
1587d6ab1dc0SHong Zhang 
1588d6ab1dc0SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1589d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1590d6ab1dc0SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]);
1591d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1592d6ab1dc0SHong Zhang 
1593d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1594d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1595d6ab1dc0SHong Zhang   /* determine row ownership */
1596d6ab1dc0SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1597d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1598d6ab1dc0SHong Zhang   merge->rowmap->n = pn;
1599d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
1600d6ab1dc0SHong Zhang   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1601d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1602d6ab1dc0SHong Zhang 
1603d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
1604d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1605d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1606d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1607d6ab1dc0SHong Zhang   len_s = merge->len_s;
1608d6ab1dc0SHong Zhang   merge->nsend = 0;
1609d6ab1dc0SHong Zhang 
1610d6ab1dc0SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1611d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1612d6ab1dc0SHong Zhang 
1613d6ab1dc0SHong Zhang   proc = 0;
1614d6ab1dc0SHong Zhang   for (i=0; i<pon; i++){
1615d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1616d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1617d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1618d6ab1dc0SHong Zhang   }
1619d6ab1dc0SHong Zhang 
1620d6ab1dc0SHong Zhang   len   = 0;  /* max length of buf_si[] */
1621d6ab1dc0SHong Zhang   owners_co[0] = 0;
1622d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++){
1623d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1624d6ab1dc0SHong Zhang     if (len_si[proc]){
1625d6ab1dc0SHong Zhang       merge->nsend++;
1626d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1627d6ab1dc0SHong Zhang       len += len_si[proc];
1628d6ab1dc0SHong Zhang     }
1629d6ab1dc0SHong Zhang   }
1630d6ab1dc0SHong Zhang 
1631d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
1632d6ab1dc0SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1633d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1634d6ab1dc0SHong Zhang 
1635d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1636d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1637d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1638d6ab1dc0SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1639d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++){
1640d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1641d6ab1dc0SHong Zhang     i = owners_co[proc];
1642d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1643d6ab1dc0SHong Zhang     k++;
1644d6ab1dc0SHong Zhang   }
1645d6ab1dc0SHong Zhang 
1646d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1647d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1648d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++){
1649c280ed6aSJed Brown     PetscMPIInt icompleted;
1650c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1651d6ab1dc0SHong Zhang   }
1652d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1653d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1654d6ab1dc0SHong Zhang 
1655d6ab1dc0SHong Zhang   /* send and recv coi */
1656d6ab1dc0SHong Zhang   /*-------------------*/
1657d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1658d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1659d6ab1dc0SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1660d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1661d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++){
1662d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1663d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1664d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1665d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1666d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1667d6ab1dc0SHong Zhang     */
1668d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1669d6ab1dc0SHong Zhang     nrows = len_si[proc]/2 - 1;
1670d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1671d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1672d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1673d6ab1dc0SHong Zhang     nrows = 0;
1674d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
1675d6ab1dc0SHong Zhang       nzi = coi[i+1] - coi[i];
1676d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1677d6ab1dc0SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
1678d6ab1dc0SHong Zhang       nrows++;
1679d6ab1dc0SHong Zhang     }
1680d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1681d6ab1dc0SHong Zhang     k++;
1682d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1683d6ab1dc0SHong Zhang   }
1684d6ab1dc0SHong Zhang   i = merge->nrecv;
1685d6ab1dc0SHong Zhang   while (i--) {
1686c280ed6aSJed Brown     PetscMPIInt icompleted;
1687c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1688d6ab1dc0SHong Zhang   }
1689d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1690d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1691d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1692d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1693d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1694d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1695d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1696d6ab1dc0SHong Zhang 
1697d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1698d6ab1dc0SHong Zhang   /*------------------------------------------*/
1699d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1700d6ab1dc0SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1701d6ab1dc0SHong Zhang   bi[0] = 0;
1702d6ab1dc0SHong Zhang 
1703d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1704d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1705d6ab1dc0SHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1706d6ab1dc0SHong Zhang   current_space = free_space;
1707d6ab1dc0SHong Zhang 
1708d6ab1dc0SHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1709d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++){
1710d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1711d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1712d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1713d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
1714d6ab1dc0SHong Zhang   }
1715d6ab1dc0SHong Zhang 
1716d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1717d6ab1dc0SHong Zhang   rmax = 0;
1718d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1719d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1720d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1721d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1722d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++){
1723d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1724d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1725d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1726d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1727d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1728d6ab1dc0SHong Zhang     }
1729d6ab1dc0SHong Zhang 
1730d6ab1dc0SHong Zhang     /* add received col data into lnk */
1731d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
1732d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1733d6ab1dc0SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
1734d6ab1dc0SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
1735d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1736d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1737d6ab1dc0SHong Zhang       }
1738d6ab1dc0SHong Zhang     }
1739d6ab1dc0SHong Zhang     nnz = lnk[0];
1740d6ab1dc0SHong Zhang 
1741d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1742d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1743d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1744d6ab1dc0SHong Zhang       nspacedouble++;
1745d6ab1dc0SHong Zhang     }
1746d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1747d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1748d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1749d6ab1dc0SHong Zhang     current_space->array           += nnz;
1750d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1751d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
1752d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1753d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1754d6ab1dc0SHong Zhang   }
1755d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1756d6ab1dc0SHong Zhang 
1757d6ab1dc0SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1758d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1759d6ab1dc0SHong Zhang   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]);
1760d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1761d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1762d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1763d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1764d6ab1dc0SHong Zhang 
1765d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1766d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
1767d6ab1dc0SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1768d6ab1dc0SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1769d6ab1dc0SHong Zhang 
1770d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1771d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1772d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1773d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1774d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1775d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1776d6ab1dc0SHong Zhang   for (i=0; i<pn; i++){
1777d6ab1dc0SHong Zhang     row = i + rstart;
1778d6ab1dc0SHong Zhang     nnz = bi[i+1] - bi[i];
1779d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1780d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1781d6ab1dc0SHong Zhang   }
1782d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1783d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1784d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1785d6ab1dc0SHong Zhang 
1786d6ab1dc0SHong Zhang   merge->bi            = bi;
1787d6ab1dc0SHong Zhang   merge->bj            = bj;
1788d6ab1dc0SHong Zhang   merge->coi           = coi;
1789d6ab1dc0SHong Zhang   merge->coj           = coj;
1790d6ab1dc0SHong Zhang   merge->buf_ri        = buf_ri;
1791d6ab1dc0SHong Zhang   merge->buf_rj        = buf_rj;
1792d6ab1dc0SHong Zhang   merge->owners_co     = owners_co;
1793d6ab1dc0SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
1794d6ab1dc0SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
1795d6ab1dc0SHong Zhang 
1796d6ab1dc0SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1797d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1798d6ab1dc0SHong Zhang 
1799d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1800d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
1801d6ab1dc0SHong Zhang   c->ptap        = ptap;
1802d6ab1dc0SHong Zhang   ptap->api      = PETSC_NULL;
1803d6ab1dc0SHong Zhang   ptap->apj      = PETSC_NULL;
1804d6ab1dc0SHong Zhang   ptap->merge    = merge;
1805d6ab1dc0SHong Zhang   ptap->rmax     = rmax;
1806e745005bSHong Zhang   ptap->apa      = PETSC_NULL;
1807f96d37fcSHong Zhang 
1808f96d37fcSHong Zhang   *C = Cmpi;
1809f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1810f96d37fcSHong Zhang   if (bi[pn] != 0) {
1811f96d37fcSHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
18121c7d5954SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1813f96d37fcSHong Zhang   } else {
1814f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1815f96d37fcSHong Zhang   }
1816f96d37fcSHong Zhang #endif
1817f96d37fcSHong Zhang   PetscFunctionReturn(0);
1818f96d37fcSHong Zhang }
1819