xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 5e5acdf2c7d3246fd94da6d18ffcb8751186fbe1)
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>
11af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
121634b032SHong Zhang 
13*5e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
14*5e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
15*5e5acdf2Sstefano_zampini #endif
16*5e5acdf2Sstefano_zampini 
172499ec78SHong Zhang #undef __FUNCT__
182499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
19150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
202499ec78SHong Zhang {
212499ec78SHong Zhang   PetscErrorCode ierr;
22*5e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
23*5e5acdf2Sstefano_zampini   const char     *algTypes[3] = {"scalable","nonscalable","hypre"};
24*5e5acdf2Sstefano_zampini   PetscInt       nalg = 3;
25*5e5acdf2Sstefano_zampini #else
260fc8cf34SHong Zhang   const char     *algTypes[2] = {"scalable","nonscalable"};
27*5e5acdf2Sstefano_zampini   PetscInt       nalg = 2;
28*5e5acdf2Sstefano_zampini #endif
290d3441aeSHong Zhang   PetscInt       alg = 1; /* set default algorithm */
30e55be12dSHong Zhang   MPI_Comm       comm;
312499ec78SHong Zhang 
322499ec78SHong Zhang   PetscFunctionBegin;
332499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
34e55be12dSHong Zhang     ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
356c4ed002SBarry Smith     if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
36e55be12dSHong Zhang 
370fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
38*5e5acdf2Sstefano_zampini     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,NULL);CHKERRQ(ierr);
390fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
400fc8cf34SHong Zhang 
413ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
420fc8cf34SHong Zhang     switch (alg) {
430fc8cf34SHong Zhang     case 1:
440fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
450fc8cf34SHong Zhang       break;
46*5e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
47*5e5acdf2Sstefano_zampini     case 2:
48*5e5acdf2Sstefano_zampini       ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
49*5e5acdf2Sstefano_zampini       break;
50*5e5acdf2Sstefano_zampini #endif
510fc8cf34SHong Zhang     default:
52b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
530fc8cf34SHong Zhang       break;
540fc8cf34SHong Zhang     }
553ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
562499ec78SHong Zhang   }
573ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
58598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
593ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
602499ec78SHong Zhang   PetscFunctionReturn(0);
612499ec78SHong Zhang }
622499ec78SHong Zhang 
63f33d1a9aSHong Zhang #undef __FUNCT__
64a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
65a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
66598bc09dSHong Zhang {
67598bc09dSHong Zhang   PetscErrorCode ierr;
68598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
69598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
70598bc09dSHong Zhang 
71598bc09dSHong Zhang   PetscFunctionBegin;
72b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
73598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
74a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
75a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
76ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
77a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
78a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
79598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
80598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
81598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
82598bc09dSHong Zhang   PetscFunctionReturn(0);
83598bc09dSHong Zhang }
84598bc09dSHong Zhang 
852499ec78SHong Zhang #undef __FUNCT__
86a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
87a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
884ae0847bSHong Zhang {
894ae0847bSHong Zhang   PetscErrorCode ierr;
904ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
914ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
924ae0847bSHong Zhang 
934ae0847bSHong Zhang   PetscFunctionBegin;
944ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
9526fbe8dcSKarl Rupp 
964ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
974ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
984ae0847bSHong Zhang   PetscFunctionReturn(0);
994ae0847bSHong Zhang }
1004ae0847bSHong Zhang 
1014ae0847bSHong Zhang #undef __FUNCT__
102f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
103f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
104598bc09dSHong Zhang {
105598bc09dSHong Zhang   PetscErrorCode ierr;
1064ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
107598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
108598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
109c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
110598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
1119ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
112c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
113598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
114c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
11529825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
116598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1179467f45fSHong Zhang   MPI_Comm       comm;
1189467f45fSHong Zhang   PetscMPIInt    size;
119598bc09dSHong Zhang 
120598bc09dSHong Zhang   PetscFunctionBegin;
1219467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1229467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1239467f45fSHong Zhang 
124a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
125598bc09dSHong Zhang   /*-----------------------------------------------------*/
126a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
127b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
128a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
129598bc09dSHong Zhang 
130598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
131598bc09dSHong Zhang   /*----------------------------------------------------------*/
132598bc09dSHong Zhang   /* get data from symbolic products */
133a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
134c12557a1SHong Zhang   p_oth = NULL;
1359467f45fSHong Zhang   if (size >1) {
1369467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1379467f45fSHong Zhang   }
138598bc09dSHong Zhang 
139598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
140598bc09dSHong Zhang   apa = ptap->apa;
141598bc09dSHong Zhang 
14229825010SHong Zhang   api = ptap->api;
14329825010SHong Zhang   apj = ptap->apj;
144598bc09dSHong Zhang   for (i=0; i<cm; i++) {
145c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
146e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
147598bc09dSHong Zhang 
148598bc09dSHong Zhang     /* set values in C */
149598bc09dSHong Zhang     apJ  = apj + api[i];
150598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
151598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
152598bc09dSHong Zhang 
153598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
154598bc09dSHong Zhang     ca = coa + co->i[i];
155598bc09dSHong Zhang     k  = 0;
156598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
157598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
158598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
159598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
160598bc09dSHong Zhang       k++;
161598bc09dSHong Zhang     }
162598bc09dSHong Zhang 
163598bc09dSHong Zhang     /* diagonal part of C */
164598bc09dSHong Zhang     ca = cda + cd->i[i];
165598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
166598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
167598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
168598bc09dSHong Zhang       k++;
169598bc09dSHong Zhang     }
170598bc09dSHong Zhang 
171598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
172598bc09dSHong Zhang     ca = coa + co->i[i];
173598bc09dSHong Zhang     for (; k0<conz; k0++) {
174598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
175598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
176598bc09dSHong Zhang       k++;
177598bc09dSHong Zhang     }
178598bc09dSHong Zhang   }
179598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
181598bc09dSHong Zhang   PetscFunctionReturn(0);
182598bc09dSHong Zhang }
183598bc09dSHong Zhang 
184598bc09dSHong Zhang #undef __FUNCT__
1850fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
1860fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
187598bc09dSHong Zhang {
188598bc09dSHong Zhang   PetscErrorCode     ierr;
189ce94432eSBarry Smith   MPI_Comm           comm;
1909467f45fSHong Zhang   PetscMPIInt        size;
191bfcd1a73SHong Zhang   Mat                Cmpi;
192598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1930298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1944ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
195bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1964ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
1974ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
198d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
199e55be12dSHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n,Crmax;
200598bc09dSHong Zhang   PetscBT            lnkbt;
201598bc09dSHong Zhang   PetscScalar        *apa;
202bfcd1a73SHong Zhang   PetscReal          afill;
203e55be12dSHong Zhang   PetscTable         ta;
204598bc09dSHong Zhang 
205598bc09dSHong Zhang   PetscFunctionBegin;
206ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
2079467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2089467f45fSHong Zhang 
209598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
210b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
211598bc09dSHong Zhang 
212598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
213b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
21419f0e57aSHong Zhang 
215598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
216a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
217598bc09dSHong Zhang 
218a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
219598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2209467f45fSHong Zhang   if (size > 1) {
2219467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
222598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
22320e3dc75SHong Zhang   } else {
22420e3dc75SHong Zhang     p_oth = NULL;
22520e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2269467f45fSHong Zhang   }
227598bc09dSHong Zhang 
228598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
229598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
230854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
231a1a4e84aSHong Zhang   ptap->api = api;
232598bc09dSHong Zhang   api[0]    = 0;
233598bc09dSHong Zhang 
234598bc09dSHong Zhang   /* create and initialize a linked list */
235c373ccc6SHong Zhang   ierr = PetscTableCreate(pN,pN,&ta);CHKERRQ(ierr);
2364b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
2374b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
238e55be12dSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
239e55be12dSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
240e55be12dSHong Zhang 
241e55be12dSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
242598bc09dSHong Zhang 
243bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
244f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
245598bc09dSHong Zhang   current_space = free_space;
246598bc09dSHong Zhang 
247598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
248598bc09dSHong Zhang   for (i=0; i<am; i++) {
249598bc09dSHong Zhang     /* diagonal portion of A */
250598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
251598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
252598bc09dSHong Zhang       row  = *adj++;
253598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
254598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
255598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2561fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
257598bc09dSHong Zhang     }
258598bc09dSHong Zhang     /* off-diagonal portion of A */
259598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
260598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
261598bc09dSHong Zhang       row  = *aoj++;
262598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
263598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2641fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
265598bc09dSHong Zhang     }
266598bc09dSHong Zhang 
267d14fa243SHong Zhang     apnz     = lnk[0];
268598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
269598bc09dSHong Zhang 
270598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
271598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
272f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
273598bc09dSHong Zhang       nspacedouble++;
274598bc09dSHong Zhang     }
275598bc09dSHong Zhang 
276598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2771fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
278598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2792205254eSKarl Rupp 
280598bc09dSHong Zhang     current_space->array           += apnz;
281598bc09dSHong Zhang     current_space->local_used      += apnz;
282598bc09dSHong Zhang     current_space->local_remaining -= apnz;
283598bc09dSHong Zhang   }
284598bc09dSHong Zhang 
285598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
286598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
287854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
288a1a4e84aSHong Zhang   apj  = ptap->apj;
289a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
290598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
291598bc09dSHong Zhang 
292f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
2931795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
2942205254eSKarl Rupp 
295f84c536bSHong Zhang   ptap->apa = apa;
296f84c536bSHong Zhang 
297bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
298598bc09dSHong Zhang   /*----------------------------------------------------*/
299bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
300bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
30133d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
302a2f3521dSMark F. Adams 
303bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
304bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
305598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
306598bc09dSHong Zhang   for (i=0; i<am; i++) {
307598bc09dSHong Zhang     row  = i + rstart;
308598bc09dSHong Zhang     apnz = api[i+1] - api[i];
309bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
310598bc09dSHong Zhang     apj += apnz;
311598bc09dSHong Zhang   }
312bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
314598bc09dSHong Zhang 
315bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
316bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
317f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
318bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
319bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
320598bc09dSHong Zhang 
321bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
322bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
323598bc09dSHong Zhang   c->ptap = ptap;
324598bc09dSHong Zhang 
325bfcd1a73SHong Zhang   *C = Cmpi;
326bfcd1a73SHong Zhang 
327bfcd1a73SHong Zhang   /* set MatInfo */
328118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
329bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
330bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
331bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
332bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
333bfcd1a73SHong Zhang 
334bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
335bfcd1a73SHong Zhang   if (api[am]) {
33657622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
33757622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
338bfcd1a73SHong Zhang   } else {
339bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
340bfcd1a73SHong Zhang   }
341bfcd1a73SHong Zhang #endif
342598bc09dSHong Zhang   PetscFunctionReturn(0);
343598bc09dSHong Zhang }
344598bc09dSHong Zhang 
3459123193aSHong Zhang #undef __FUNCT__
3469123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
347150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3489123193aSHong Zhang {
3499123193aSHong Zhang   PetscErrorCode ierr;
3509123193aSHong Zhang 
3519123193aSHong Zhang   PetscFunctionBegin;
3529123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3533ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3549123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3553ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3569123193aSHong Zhang   }
3573ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3589123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3593ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3609123193aSHong Zhang   PetscFunctionReturn(0);
3619123193aSHong Zhang }
3629123193aSHong Zhang 
3634324174eSBarry Smith typedef struct {
3644324174eSBarry Smith   Mat         workB;
3654324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3664324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3674324174eSBarry Smith } MPIAIJ_MPIDense;
3684324174eSBarry Smith 
3697af9e4fdSHong Zhang #undef __FUNCT__
37096b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
37196b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3724324174eSBarry Smith {
3734324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3744324174eSBarry Smith   PetscErrorCode  ierr;
3754324174eSBarry Smith 
3764324174eSBarry Smith   PetscFunctionBegin;
3776bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3784324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3794324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3804324174eSBarry Smith   PetscFunctionReturn(0);
3814324174eSBarry Smith }
3824324174eSBarry Smith 
3839123193aSHong Zhang #undef __FUNCT__
384fd4e9aacSBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense"
385fd4e9aacSBarry Smith /*
386fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
387fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
388fd4e9aacSBarry Smith 
389fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
390fd4e9aacSBarry Smith */
391fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
392fd4e9aacSBarry Smith {
393fd4e9aacSBarry Smith   PetscErrorCode         ierr;
394fd4e9aacSBarry Smith   PetscBool              flg;
395fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
396fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
397fd4e9aacSBarry Smith   PetscContainer         container;
398fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
399fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
400fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
401fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
402fd4e9aacSBarry Smith 
403fd4e9aacSBarry Smith   PetscFunctionBegin;
404fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
405fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
406fd4e9aacSBarry Smith 
407fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
408fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
409fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
410fd4e9aacSBarry Smith 
411fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
412fd4e9aacSBarry Smith 
413fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
414fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
415fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
416fd4e9aacSBarry Smith   /* Create work arrays needed */
417fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
418fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
419fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
420fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
421fd4e9aacSBarry Smith 
422fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
423fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
424fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
425fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
426fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
427fd4e9aacSBarry Smith 
428fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
429fd4e9aacSBarry Smith   PetscFunctionReturn(0);
430fd4e9aacSBarry Smith }
431fd4e9aacSBarry Smith 
432fd4e9aacSBarry Smith #undef __FUNCT__
4339123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4349123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4359123193aSHong Zhang {
4369123193aSHong Zhang   PetscErrorCode         ierr;
437f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
438d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
439bf0cc555SLisandro Dalcin   PetscContainer         container;
4404324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4414324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4424324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4434324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
444d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4459123193aSHong Zhang 
4469123193aSHong Zhang   PetscFunctionBegin;
447ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
448d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
44933d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
450cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4510298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
452cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
453cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4542205254eSKarl Rupp 
455d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
456f32d5d43SBarry Smith 
457b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
458f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4590298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4604324174eSBarry Smith   /* Create work arrays needed */
461dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
462dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
463dcca6d9dSJed Brown                       from->n,&contents->rwaits,
464dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4654324174eSBarry Smith 
466ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
467bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
46896b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
469bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
470bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4719123193aSHong Zhang   PetscFunctionReturn(0);
4729123193aSHong Zhang }
4739123193aSHong Zhang 
4747af9e4fdSHong Zhang #undef __FUNCT__
4757af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
476f32d5d43SBarry Smith /*
4772636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4782636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
479f32d5d43SBarry Smith */
4804324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
481f32d5d43SBarry Smith {
482f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
483f32d5d43SBarry Smith   PetscErrorCode         ierr;
484f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
485f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
486f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
487f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
488f32d5d43SBarry Smith   PetscInt               i,j,k;
489aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
490aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
491f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
492ce94432eSBarry Smith   MPI_Comm               comm;
493d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
494f32d5d43SBarry Smith   MPI_Status             status;
4954324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
496bf0cc555SLisandro Dalcin   PetscContainer         container;
4974324174eSBarry Smith   Mat                    workB;
498f32d5d43SBarry Smith 
499f32d5d43SBarry Smith   PetscFunctionBegin;
500ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
501bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
50229825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
503bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
5044324174eSBarry Smith 
5054324174eSBarry Smith   workB = *outworkB = contents->workB;
506cf1a0672SHong 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);
507f32d5d43SBarry Smith   sindices = to->indices;
508f32d5d43SBarry Smith   sstarts  = to->starts;
509f32d5d43SBarry Smith   sprocs   = to->procs;
5104324174eSBarry Smith   swaits   = contents->swaits;
5114324174eSBarry Smith   svalues  = contents->svalues;
512f32d5d43SBarry Smith 
513f32d5d43SBarry Smith   rindices = from->indices;
514f32d5d43SBarry Smith   rstarts  = from->starts;
515f32d5d43SBarry Smith   rprocs   = from->procs;
5164324174eSBarry Smith   rwaits   = contents->rwaits;
5174324174eSBarry Smith   rvalues  = contents->rvalues;
518f32d5d43SBarry Smith 
5198c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5208c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
521f32d5d43SBarry Smith 
522f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
523f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
524f32d5d43SBarry Smith   }
525f32d5d43SBarry Smith 
526f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
527f32d5d43SBarry Smith     /* pack a message at a time */
528f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
529f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
530f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5312636ff26SBarry Smith       }
5322636ff26SBarry Smith     }
533f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
534f32d5d43SBarry Smith   }
5352636ff26SBarry Smith 
536f32d5d43SBarry Smith   nrecvs = from->n;
537f32d5d43SBarry Smith   while (nrecvs) {
538f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
539f32d5d43SBarry Smith     nrecvs--;
540f32d5d43SBarry Smith     /* unpack a message at a time */
541f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
542f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
543f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5442636ff26SBarry Smith       }
5452636ff26SBarry Smith     }
546f32d5d43SBarry Smith   }
547cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
548f32d5d43SBarry Smith 
5498c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5508c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
551f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
552f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
553f32d5d43SBarry Smith   PetscFunctionReturn(0);
554f32d5d43SBarry Smith }
555f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
556f32d5d43SBarry Smith 
5579123193aSHong Zhang #undef __FUNCT__
5589123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5599123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5609123193aSHong Zhang {
5619123193aSHong Zhang   PetscErrorCode ierr;
562f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
563f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
564f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
565f32d5d43SBarry Smith   Mat            workB;
5669123193aSHong Zhang 
5679123193aSHong Zhang   PetscFunctionBegin;
568f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
569f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
570f32d5d43SBarry Smith 
571f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5724324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
573f32d5d43SBarry Smith 
574f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
575f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5769123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5779123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5789123193aSHong Zhang   PetscFunctionReturn(0);
5799123193aSHong Zhang }
580cf1a0672SHong Zhang 
5811634b032SHong Zhang #undef __FUNCT__
582f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
583f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5841634b032SHong Zhang {
5851634b032SHong Zhang   PetscErrorCode ierr;
586cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
587cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
588cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
589cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
590cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
591cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
592cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
59329825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
59429825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
595cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
59629825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
597cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
59829825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
59929825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
600a7c7454dSHong Zhang   MPI_Comm       comm;
601a7c7454dSHong Zhang   PetscMPIInt    size;
6021634b032SHong Zhang 
6031634b032SHong Zhang   PetscFunctionBegin;
604a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
605a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
606a7c7454dSHong Zhang 
607cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
608cf1a0672SHong Zhang   /*-----------------------------------------------------*/
609cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
610b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
611cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
612cf1a0672SHong Zhang 
613cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
614cf1a0672SHong Zhang   /*----------------------------------------------------------*/
615cf1a0672SHong Zhang   /* get data from symbolic products */
616cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
617cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
618a7c7454dSHong Zhang   if (size >1) {
619a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
620cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
62120e3dc75SHong Zhang   } else {
62220e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
623a7c7454dSHong Zhang   }
624cf1a0672SHong Zhang 
62529825010SHong Zhang   api = ptap->api;
62629825010SHong Zhang   apj = ptap->apj;
627cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
62829825010SHong Zhang     apJ = apj + api[i];
62929825010SHong Zhang 
630cf1a0672SHong Zhang     /* diagonal portion of A */
631cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
632cf1a0672SHong Zhang     adj = ad->j + adi[i];
633cf1a0672SHong Zhang     ada = ad->a + adi[i];
634cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
635cf1a0672SHong Zhang       row = adj[j];
636cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
637cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
638cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
63929825010SHong Zhang       /* perform sparse axpy */
640cf1a0672SHong Zhang       valtmp = ada[j];
64129825010SHong Zhang       nextp  = 0;
64229825010SHong Zhang       for (k=0; nextp<pnz; k++) {
64329825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
64429825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
64529825010SHong Zhang         }
6461634b032SHong Zhang       }
647cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
648cf1a0672SHong Zhang     }
6491634b032SHong Zhang 
650cf1a0672SHong Zhang     /* off-diagonal portion of A */
651cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
652cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
653cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
654cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
655cf1a0672SHong Zhang       row = aoj[j];
656cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
657cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
658cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
65929825010SHong Zhang       /* perform sparse axpy */
660cf1a0672SHong Zhang       valtmp = aoa[j];
66129825010SHong Zhang       nextp  = 0;
66229825010SHong Zhang       for (k=0; nextp<pnz; k++) {
66329825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
66429825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
66529825010SHong Zhang         }
666cf1a0672SHong Zhang       }
667cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
668cf1a0672SHong Zhang     }
6691634b032SHong Zhang 
670cf1a0672SHong Zhang     /* set values in C */
671cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
672cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6731634b032SHong Zhang 
674cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
675cf1a0672SHong Zhang     ca = coa + co->i[i];
676cf1a0672SHong Zhang     k  = 0;
677cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
678cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
67929825010SHong Zhang       ca[k0]        = apa_sparse[k];
68029825010SHong Zhang       apa_sparse[k] = 0.0;
681cf1a0672SHong Zhang       k++;
682cf1a0672SHong Zhang     }
6831634b032SHong Zhang 
684cf1a0672SHong Zhang     /* diagonal part of C */
685cf1a0672SHong Zhang     ca = cda + cd->i[i];
686cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
68729825010SHong Zhang       ca[k1]        = apa_sparse[k];
68829825010SHong Zhang       apa_sparse[k] = 0.0;
689cf1a0672SHong Zhang       k++;
690cf1a0672SHong Zhang     }
691cf1a0672SHong Zhang 
692cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
693cf1a0672SHong Zhang     ca = coa + co->i[i];
694cf1a0672SHong Zhang     for (; k0<conz; k0++) {
69529825010SHong Zhang       ca[k0]        = apa_sparse[k];
69629825010SHong Zhang       apa_sparse[k] = 0.0;
697cf1a0672SHong Zhang       k++;
698cf1a0672SHong Zhang     }
699cf1a0672SHong Zhang   }
700cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
701cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
702cf1a0672SHong Zhang   PetscFunctionReturn(0);
703cf1a0672SHong Zhang }
704cf1a0672SHong Zhang 
7050fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
706cf1a0672SHong Zhang #undef __FUNCT__
707b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
708b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
709cf1a0672SHong Zhang {
710cf1a0672SHong Zhang   PetscErrorCode     ierr;
711ce94432eSBarry Smith   MPI_Comm           comm;
712a7c7454dSHong Zhang   PetscMPIInt        size;
713cf1a0672SHong Zhang   Mat                Cmpi;
714cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
7150298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
716cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
717cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
718cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
719cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
7200ca7d551SHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max;
7210ca7d551SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
722cf1a0672SHong Zhang   PetscReal          afill;
723cf1a0672SHong Zhang   PetscScalar        *apa;
724eca6b86aSHong Zhang   PetscTable         ta;
725cf1a0672SHong Zhang 
726cf1a0672SHong Zhang   PetscFunctionBegin;
727ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
728a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
729a7c7454dSHong Zhang 
730cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
731b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
732cf1a0672SHong Zhang 
733cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
734b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
73519f0e57aSHong Zhang 
736cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
737cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
738cf1a0672SHong Zhang 
739cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
740cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
741a7c7454dSHong Zhang   if (size > 1) {
742a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
74320e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
744968382fdSHong Zhang   } else {
74520e3dc75SHong Zhang     p_oth  = NULL;
74620e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
747a7c7454dSHong Zhang   }
748cf1a0672SHong Zhang 
749cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
750cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
751854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
752cf1a0672SHong Zhang   ptap->api = api;
753cf1a0672SHong Zhang   api[0]    = 0;
754cf1a0672SHong Zhang 
755cf1a0672SHong Zhang   /* create and initialize a linked list */
7560ca7d551SHong Zhang   apnz_max = 6*(p_loc->rmax + (PetscInt)(1.e-2*pN)); /* expected apnz_max */
7570ca7d551SHong Zhang   if (apnz_max > pN) apnz_max = pN;
7580ca7d551SHong Zhang   ierr = PetscTableCreate(apnz_max,pN,&ta);CHKERRQ(ierr);
7590ca7d551SHong Zhang 
7600ca7d551SHong Zhang   /* Calculate apnz_max */
7610ca7d551SHong Zhang   apnz_max = 0;
7620ca7d551SHong Zhang   for (i=0; i<am; i++) {
7630ca7d551SHong Zhang     ierr = PetscTableRemoveAll(ta);CHKERRQ(ierr);
7640ca7d551SHong Zhang     /* diagonal portion of A */
7650ca7d551SHong Zhang     nzi  = adi[i+1] - adi[i];
7660ca7d551SHong Zhang     Jptr = adj+adi[i];  /* cols of A_diag */
7670ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_loc,nzi,Jptr,ta);
7680ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7690ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7700ca7d551SHong Zhang 
7710ca7d551SHong Zhang     /*  off-diagonal portion of A */
7720ca7d551SHong Zhang     nzi = aoi[i+1] - aoi[i];
7730ca7d551SHong Zhang     Jptr = aoj+aoi[i];  /* cols of A_off */
7740ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_oth,nzi,Jptr,ta);
7750ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7760ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7770ca7d551SHong Zhang   }
778eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
779eca6b86aSHong Zhang 
7800ca7d551SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(apnz_max,&lnk);CHKERRQ(ierr);
781cf1a0672SHong Zhang 
782cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
783f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
784cf1a0672SHong Zhang   current_space = free_space;
785cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
786cf1a0672SHong Zhang   for (i=0; i<am; i++) {
787cf1a0672SHong Zhang     /* diagonal portion of A */
788cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
789cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
790cf1a0672SHong Zhang       row  = *adj++;
791cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
792cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
793cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
794f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
795cf1a0672SHong Zhang     }
796cf1a0672SHong Zhang     /* off-diagonal portion of A */
797cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
798cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
799cf1a0672SHong Zhang       row  = *aoj++;
800cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
801cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
802f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
803cf1a0672SHong Zhang     }
804cf1a0672SHong Zhang 
805f84c536bSHong Zhang     apnz     = *lnk;
806cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
807cf1a0672SHong Zhang 
808cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
809cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
810f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
811cf1a0672SHong Zhang       nspacedouble++;
812cf1a0672SHong Zhang     }
813cf1a0672SHong Zhang 
814cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
815f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
816cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
8172205254eSKarl Rupp 
818cf1a0672SHong Zhang     current_space->array           += apnz;
819cf1a0672SHong Zhang     current_space->local_used      += apnz;
820cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
821cf1a0672SHong Zhang   }
822cf1a0672SHong Zhang 
823cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
824cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
825854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
826cf1a0672SHong Zhang   apj  = ptap->apj;
827cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
828f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
829cf1a0672SHong Zhang 
830cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
831cf1a0672SHong Zhang   /*----------------------------------------------------*/
832cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
833cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
83433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
835cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
836cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
837cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
838cf1a0672SHong Zhang 
83929825010SHong Zhang   /* malloc apa for assembly Cmpi */
8401795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8412205254eSKarl Rupp 
84229825010SHong Zhang   ptap->apa = apa;
843cf1a0672SHong Zhang   for (i=0; i<am; i++) {
844cf1a0672SHong Zhang     row  = i + rstart;
845cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
846cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
847cf1a0672SHong Zhang     apj += apnz;
848cf1a0672SHong Zhang   }
849cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
850cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
851cf1a0672SHong Zhang 
852cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
853cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
854f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
855cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
856cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
857cf1a0672SHong Zhang 
858cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
859cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
860cf1a0672SHong Zhang   c->ptap = ptap;
861cf1a0672SHong Zhang 
862cf1a0672SHong Zhang   *C = Cmpi;
863cf1a0672SHong Zhang 
864cf1a0672SHong Zhang   /* set MatInfo */
865118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
866cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
867cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
868cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
869cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
870cf1a0672SHong Zhang 
871cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
872cf1a0672SHong Zhang   if (api[am]) {
87357622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
87457622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
875cf1a0672SHong Zhang   } else {
876cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
877cf1a0672SHong Zhang   }
878cf1a0672SHong Zhang #endif
8791634b032SHong Zhang   PetscFunctionReturn(0);
8801634b032SHong Zhang }
881f96d37fcSHong Zhang 
882f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
883f96d37fcSHong Zhang #undef __FUNCT__
884f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
885c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
886f96d37fcSHong Zhang {
887f96d37fcSHong Zhang   PetscErrorCode ierr;
888c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
889c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
890f96d37fcSHong Zhang 
891f96d37fcSHong Zhang   PetscFunctionBegin;
892c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
893c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
894c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
895c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
896c216dbf3SHong Zhang 
897c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
898c216dbf3SHong Zhang     switch (alg) {
899c216dbf3SHong Zhang     case 1:
9002bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
901c216dbf3SHong Zhang       break;
902c216dbf3SHong Zhang     case 2:
903275476c6SMatthew G. Knepley     {
904ab1b013aSHong Zhang       Mat         Pt;
905ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
906ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
907ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
908ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
909ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
910ab1b013aSHong Zhang       ptap     = c->ptap;
911ab1b013aSHong Zhang       ptap->Pt = Pt;
912c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
913c216dbf3SHong Zhang       PetscFunctionReturn(0);
914275476c6SMatthew G. Knepley     }
915c216dbf3SHong Zhang       break;
916c216dbf3SHong Zhang     default:
9176da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
918c216dbf3SHong Zhang       break;
919c216dbf3SHong Zhang     }
920c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
921c216dbf3SHong Zhang   }
922c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
923c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
924c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
925ab1b013aSHong Zhang   PetscFunctionReturn(0);
926ab1b013aSHong Zhang }
927ab1b013aSHong Zhang 
928c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
929c216dbf3SHong Zhang #undef __FUNCT__
930c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
931c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
932c216dbf3SHong Zhang {
933c216dbf3SHong Zhang   PetscErrorCode ierr;
9342bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9352bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9362bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
937c216dbf3SHong Zhang 
938c216dbf3SHong Zhang   PetscFunctionBegin;
939c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
940c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
941f96d37fcSHong Zhang   PetscFunctionReturn(0);
942f96d37fcSHong Zhang }
943f96d37fcSHong Zhang 
9446da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
945f96d37fcSHong Zhang #undef __FUNCT__
9466da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9476da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
948f96d37fcSHong Zhang {
949c5af039cSHong Zhang   PetscErrorCode      ierr;
950c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
951497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
952c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
953c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
954d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
955497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
956e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
957c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
958ce94432eSBarry Smith   MPI_Comm            comm;
959c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
960c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
961c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
962c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
963c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
964c5af039cSHong Zhang   MPI_Status          *status;
965c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
966d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
967a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
968c5af039cSHong Zhang   Mat                 A_loc;
969c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
970f96d37fcSHong Zhang 
971f96d37fcSHong Zhang   PetscFunctionBegin;
972ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
973c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
974c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
975c5af039cSHong Zhang 
976c5af039cSHong Zhang   ptap  = c->ptap;
977c5af039cSHong Zhang   merge = ptap->merge;
978c5af039cSHong Zhang 
979c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
980c5af039cSHong Zhang   /*--------------------------------------------------------------*/
981c5af039cSHong Zhang   /* get data from symbolic products */
982c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
983854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
984c5af039cSHong Zhang 
985c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
986c5af039cSHong Zhang   owners = merge->rowmap->range;
987854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
988c5af039cSHong Zhang 
989c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
990c5af039cSHong Zhang   A_loc = ptap->A_loc;
991c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
992c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
993d6ab1dc0SHong Zhang   ai    = a_loc->i;
994d6ab1dc0SHong Zhang   aj    = a_loc->j;
995c5af039cSHong Zhang 
996854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
997c5af039cSHong Zhang 
998c5af039cSHong Zhang   for (i=0; i<am; i++) {
999e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
1000d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1001d6ab1dc0SHong Zhang     adj = aj + ai[i];
1002d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1003c5af039cSHong Zhang     for (j=0; j<anz; j++) {
1004e745005bSHong Zhang       aval[adj[j]] = ada[j];
1005c5af039cSHong Zhang     }
1006c5af039cSHong Zhang 
1007c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1008c5af039cSHong Zhang     /*--------------------------------------------------------------*/
1009c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1010c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
1011c5af039cSHong Zhang     poJ = po->j + po->i[i];
1012c5af039cSHong Zhang     pA  = po->a + po->i[i];
1013c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1014c5af039cSHong Zhang       row = poJ[j];
1015c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
1016c5af039cSHong Zhang       cj  = coj + coi[row];
1017c5af039cSHong Zhang       ca  = coa + coi[row];
1018c5af039cSHong Zhang       /* perform dense axpy */
1019c5af039cSHong Zhang       valtmp = pA[j];
1020c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1021e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1022c5af039cSHong Zhang       }
1023c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1024c5af039cSHong Zhang     }
1025c5af039cSHong Zhang 
1026c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
1027c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1028c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1029c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1030c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1031c5af039cSHong Zhang       row = pdJ[j];
1032c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1033c5af039cSHong Zhang       cj  = bj + bi[row];
1034c5af039cSHong Zhang       ca  = ba + bi[row];
1035c5af039cSHong Zhang       /* perform dense axpy */
1036c5af039cSHong Zhang       valtmp = pA[j];
1037c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1038e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1039c5af039cSHong Zhang       }
1040c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1041c5af039cSHong Zhang     }
1042c5af039cSHong Zhang 
1043d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1044d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1045e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1046c5af039cSHong Zhang   }
1047c5af039cSHong Zhang 
1048c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1049c5af039cSHong Zhang   /*------------------------------------*/
1050c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1051c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1052c5af039cSHong Zhang   len_s  = merge->len_s;
1053c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1054c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1055c5af039cSHong Zhang 
1056dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1057c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1058c5af039cSHong Zhang     if (!len_s[proc]) continue;
1059c5af039cSHong Zhang     i    = merge->owners_co[proc];
1060c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1061c5af039cSHong Zhang     k++;
1062c5af039cSHong Zhang   }
1063c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1064c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1065c5af039cSHong Zhang 
1066c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1067c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1068c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1069c5af039cSHong Zhang 
1070c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1071c5af039cSHong Zhang   /*----------------------------------------------------*/
1072dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1073c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1074c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1075c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1076c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1077c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1078c5af039cSHong Zhang   }
1079c5af039cSHong Zhang 
1080c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1081c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1082c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1083c5af039cSHong Zhang     ba_i = ba + bi[i];
1084c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1085c5af039cSHong Zhang     /* add received vals into ba */
1086c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1087c5af039cSHong Zhang       /* i-th row */
1088c5af039cSHong Zhang       if (i == *nextrow[k]) {
1089c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1090c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1091c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1092c5af039cSHong Zhang         nextcj = 0;
1093c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1094c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1095c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1096c5af039cSHong Zhang           }
1097c5af039cSHong Zhang         }
1098c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1099c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1100c5af039cSHong Zhang       }
1101c5af039cSHong Zhang     }
1102c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1103c5af039cSHong Zhang   }
1104c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1105c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1106c5af039cSHong Zhang 
1107c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1108c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1109c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1110c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1111e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1112f96d37fcSHong Zhang   PetscFunctionReturn(0);
1113f96d37fcSHong Zhang }
1114f96d37fcSHong Zhang 
1115c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1116c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1117f96d37fcSHong Zhang #undef __FUNCT__
11182bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
11192bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1120f96d37fcSHong Zhang {
1121f96d37fcSHong Zhang   PetscErrorCode      ierr;
11224e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1123f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
11240298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1125f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1126f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1127f96d37fcSHong Zhang   PetscInt            nnz;
11284e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1129497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1130f96d37fcSHong Zhang   PetscBT             lnkbt;
1131ce94432eSBarry Smith   MPI_Comm            comm;
1132f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1133f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1134f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1135f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1136f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1137f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1138f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1139f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1140d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1141f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1142438d860cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1143f96d37fcSHong Zhang   PetscScalar         *vals;
11444e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1145f96d37fcSHong Zhang 
1146f96d37fcSHong Zhang   PetscFunctionBegin;
1147ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1148c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
11496c4ed002SBarry Smith   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1150c5af039cSHong Zhang 
1151f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1152f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1153f96d37fcSHong Zhang 
1154f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1155b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1156f96d37fcSHong Zhang 
1157f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1158f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11592205254eSKarl Rupp 
1160c5af039cSHong Zhang   ptap->A_loc = A_loc;
11612205254eSKarl Rupp 
11621c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1163d6ab1dc0SHong Zhang   ai    = a_loc->i;
1164d6ab1dc0SHong Zhang   aj    = a_loc->j;
1165f96d37fcSHong Zhang 
1166f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1167f96d37fcSHong Zhang   /*----------------------------------------------------*/
11684e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11694e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11704e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11714e938277SHong Zhang 
11724e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11734e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11744e938277SHong Zhang   poti = pot->i; potj = pot->j;
1175f96d37fcSHong Zhang 
1176f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11772205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1178854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1179f96d37fcSHong Zhang   coi[0] = 0;
1180f96d37fcSHong Zhang 
1181f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1182f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1183a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1184f96d37fcSHong Zhang   current_space = free_space;
118519f0e57aSHong Zhang 
1186f96d37fcSHong Zhang   /* create and initialize a linked list */
1187438d860cSHong Zhang   ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1188f96d37fcSHong Zhang 
1189f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1190f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1191f96d37fcSHong Zhang     ptJ = potj + poti[i];
1192f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1193f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1194d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1195d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1196f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1197d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1198f96d37fcSHong Zhang     }
11994e938277SHong Zhang     nnz = lnk[0];
1200f96d37fcSHong Zhang 
1201f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1202f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1203f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1204f96d37fcSHong Zhang       nspacedouble++;
1205f96d37fcSHong Zhang     }
1206f96d37fcSHong Zhang 
1207f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
12084e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
12092205254eSKarl Rupp 
1210f96d37fcSHong Zhang     current_space->array           += nnz;
1211f96d37fcSHong Zhang     current_space->local_used      += nnz;
1212f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
12132205254eSKarl Rupp 
1214f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1215f96d37fcSHong Zhang   }
1216f96d37fcSHong Zhang 
1217854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1218f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
12192205254eSKarl Rupp 
1220118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1221f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1222f96d37fcSHong Zhang 
1223f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1224f96d37fcSHong Zhang   /*----------------------------------------------*/
1225f96d37fcSHong Zhang   /* determine row ownership */
1226b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1227f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12282205254eSKarl Rupp 
1229f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1230f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12312205254eSKarl Rupp 
1232f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1233f96d37fcSHong Zhang   owners = merge->rowmap->range;
1234f96d37fcSHong Zhang 
1235f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12361795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1237785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12382205254eSKarl Rupp 
1239f96d37fcSHong Zhang   len_s        = merge->len_s;
1240f96d37fcSHong Zhang   merge->nsend = 0;
1241f96d37fcSHong Zhang 
1242854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1243f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1244f96d37fcSHong Zhang 
1245f96d37fcSHong Zhang   proc = 0;
1246f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1247f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1248f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1249f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1250f96d37fcSHong Zhang   }
1251f96d37fcSHong Zhang 
1252f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1253f96d37fcSHong Zhang   owners_co[0] = 0;
1254f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1255f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1256f96d37fcSHong Zhang     if (len_si[proc]) {
1257f96d37fcSHong Zhang       merge->nsend++;
1258f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1259f96d37fcSHong Zhang       len         += len_si[proc];
1260f96d37fcSHong Zhang     }
1261f96d37fcSHong Zhang   }
1262f96d37fcSHong Zhang 
1263f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12640298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1265f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1266f96d37fcSHong Zhang 
1267f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1268f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1269f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1270854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1271f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1272f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1273f96d37fcSHong Zhang     i    = owners_co[proc];
1274f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1275f96d37fcSHong Zhang     k++;
1276f96d37fcSHong Zhang   }
1277f96d37fcSHong Zhang 
1278f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1279785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1280f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1281c280ed6aSJed Brown     PetscMPIInt icompleted;
1282c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1283f96d37fcSHong Zhang   }
1284f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1285f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1286f96d37fcSHong Zhang 
1287f96d37fcSHong Zhang   /* send and recv coi */
1288f96d37fcSHong Zhang   /*-------------------*/
1289f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1290f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1291854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1292f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1293f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1294f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1295f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1296f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1297f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1298f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1299f96d37fcSHong Zhang     */
1300f96d37fcSHong Zhang     /*-------------------------------------------*/
1301f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1302f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1303f96d37fcSHong Zhang     buf_si[0]   = nrows;
1304f96d37fcSHong Zhang     buf_si_i[0] = 0;
1305f96d37fcSHong Zhang     nrows       = 0;
1306f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1307f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1308f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1309f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1310f96d37fcSHong Zhang       nrows++;
1311f96d37fcSHong Zhang     }
1312f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1313f96d37fcSHong Zhang     k++;
1314f96d37fcSHong Zhang     buf_si += len_si[proc];
1315f96d37fcSHong Zhang   }
1316f96d37fcSHong Zhang   i = merge->nrecv;
1317f96d37fcSHong Zhang   while (i--) {
1318c280ed6aSJed Brown     PetscMPIInt icompleted;
1319c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1320f96d37fcSHong Zhang   }
1321f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1322f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1323f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1324f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1325f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1326f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1327f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1328f96d37fcSHong Zhang 
1329f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1330f96d37fcSHong Zhang   /*------------------------------------------*/
1331f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1332854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1333f96d37fcSHong Zhang   bi[0] = 0;
1334f96d37fcSHong Zhang 
1335c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1336f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1337a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1338f96d37fcSHong Zhang   current_space = free_space;
1339f96d37fcSHong Zhang 
1340dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1341f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1342f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1343f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1344f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1345f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1346f96d37fcSHong Zhang   }
1347f96d37fcSHong Zhang 
13481c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1349f96d37fcSHong Zhang   rmax = 0;
1350f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1351f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1352f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1353f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1354f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1355f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1356d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1357d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1358f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1359d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1360f96d37fcSHong Zhang     }
13614e938277SHong Zhang 
1362f96d37fcSHong Zhang     /* add received col data into lnk */
1363f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1364f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1365f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1366f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13674e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1368f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1369f96d37fcSHong Zhang       }
1370f96d37fcSHong Zhang     }
13714e938277SHong Zhang     nnz = lnk[0];
1372f96d37fcSHong Zhang 
1373f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1374f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1375f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1376f96d37fcSHong Zhang       nspacedouble++;
1377f96d37fcSHong Zhang     }
1378f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13794e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1380f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13812205254eSKarl Rupp 
1382f96d37fcSHong Zhang     current_space->array           += nnz;
1383f96d37fcSHong Zhang     current_space->local_used      += nnz;
1384f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13852205254eSKarl Rupp 
1386f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1387f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1388f96d37fcSHong Zhang   }
1389f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1390f96d37fcSHong Zhang 
1391854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1392f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13932205254eSKarl Rupp 
1394118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1395f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1396d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13974e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13984e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1399f96d37fcSHong Zhang 
14001c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
14011c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
14021795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1403f96d37fcSHong Zhang 
1404f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1405f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
140633d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1407f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1408f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1409f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1410f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1411f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1412f96d37fcSHong Zhang     row  = i + rstart;
1413f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1414f96d37fcSHong Zhang     Jptr = bj + bi[i];
1415f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1416f96d37fcSHong Zhang   }
1417f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1418f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1419f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1420f96d37fcSHong Zhang 
1421f96d37fcSHong Zhang   merge->bi        = bi;
1422f96d37fcSHong Zhang   merge->bj        = bj;
1423f96d37fcSHong Zhang   merge->coi       = coi;
1424f96d37fcSHong Zhang   merge->coj       = coj;
1425f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1426f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1427f96d37fcSHong Zhang   merge->owners_co = owners_co;
1428f96d37fcSHong Zhang 
1429f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1430f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1431f96d37fcSHong Zhang   c->ptap     = ptap;
14320298fd71SBarry Smith   ptap->api   = NULL;
14330298fd71SBarry Smith   ptap->apj   = NULL;
1434f96d37fcSHong Zhang   ptap->merge = merge;
1435375ed354SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1436375ed354SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1437375ed354SHong Zhang 
1438375ed354SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1439375ed354SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1440375ed354SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1441d6ab1dc0SHong Zhang 
1442d6ab1dc0SHong Zhang   *C = Cmpi;
1443d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1444d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
144557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
144657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1447d6ab1dc0SHong Zhang   } else {
1448d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1449d6ab1dc0SHong Zhang   }
1450d6ab1dc0SHong Zhang #endif
1451d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1452d6ab1dc0SHong Zhang }
1453d6ab1dc0SHong Zhang 
1454d6ab1dc0SHong Zhang #undef __FUNCT__
14556da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14566da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1457d6ab1dc0SHong Zhang {
1458d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1459d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1460d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1461d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1462d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1463e745005bSHong Zhang   PetscInt            *adj;
1464e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1465e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1466d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1467ce94432eSBarry Smith   MPI_Comm            comm;
1468d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1469d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1470d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1471d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1472d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1473d6ab1dc0SHong Zhang   MPI_Status          *status;
1474d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1475d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1476a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1477d6ab1dc0SHong Zhang   Mat                 A_loc;
1478d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1479d6ab1dc0SHong Zhang 
1480d6ab1dc0SHong Zhang   PetscFunctionBegin;
1481ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1482d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1483d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1484d6ab1dc0SHong Zhang 
1485d6ab1dc0SHong Zhang   ptap  = c->ptap;
1486d6ab1dc0SHong Zhang   merge = ptap->merge;
1487d6ab1dc0SHong Zhang 
1488e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1489e745005bSHong Zhang   /*------------------------------------------*/
1490d6ab1dc0SHong Zhang   /* get data from symbolic products */
1491d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1492854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1493d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1494d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14951795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1496d6ab1dc0SHong Zhang 
1497d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1498d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1499d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1500d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1501d6ab1dc0SHong Zhang   ai    = a_loc->i;
1502d6ab1dc0SHong Zhang   aj    = a_loc->j;
1503d6ab1dc0SHong Zhang 
1504d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1505d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1506d6ab1dc0SHong Zhang     adj = aj + ai[i];
1507d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1508d6ab1dc0SHong Zhang 
1509d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1510e745005bSHong Zhang     /*-------------------------------------------------------------*/
1511d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1512d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1513d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1514d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1515d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1516d6ab1dc0SHong Zhang       row = poJ[j];
1517d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1518d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1519e745005bSHong Zhang       /* perform sparse axpy */
1520e745005bSHong Zhang       nexta  = 0;
1521d6ab1dc0SHong Zhang       valtmp = pA[j];
1522e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1523e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1524e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1525e745005bSHong Zhang           nexta++;
1526d6ab1dc0SHong Zhang         }
1527e745005bSHong Zhang       }
1528e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang     }
1530d6ab1dc0SHong Zhang 
1531d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1532d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1533d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1534d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1535d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1536d6ab1dc0SHong Zhang       row = pdJ[j];
1537d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1538d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1539e745005bSHong Zhang       /* perform sparse axpy */
1540e745005bSHong Zhang       nexta  = 0;
1541d6ab1dc0SHong Zhang       valtmp = pA[j];
1542e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1543e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1544e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1545e745005bSHong Zhang           nexta++;
1546d6ab1dc0SHong Zhang         }
1547e745005bSHong Zhang       }
1548e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1549d6ab1dc0SHong Zhang     }
1550d6ab1dc0SHong Zhang   }
1551d6ab1dc0SHong Zhang 
1552d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1553d6ab1dc0SHong Zhang   /*------------------------------------*/
1554d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1555d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1556d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1557d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1558d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1559d6ab1dc0SHong Zhang 
1560dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1561d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1562d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1563d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1564d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1565d6ab1dc0SHong Zhang     k++;
1566d6ab1dc0SHong Zhang   }
1567d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1568d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1569d6ab1dc0SHong Zhang 
1570d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1571d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1572d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1573d6ab1dc0SHong Zhang 
1574d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1575d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1576dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1577d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1578e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1579d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1580d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1581d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1582d6ab1dc0SHong Zhang   }
1583d6ab1dc0SHong Zhang 
1584d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1585d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1586d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1587d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1588d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1589d6ab1dc0SHong Zhang     /* add received vals into ba */
1590d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1591d6ab1dc0SHong Zhang       /* i-th row */
1592d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1593d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1594d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1595d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1596d6ab1dc0SHong Zhang         nextcj = 0;
1597d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1598d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1599d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1600d6ab1dc0SHong Zhang           }
1601d6ab1dc0SHong Zhang         }
1602d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1603d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1604d6ab1dc0SHong Zhang       }
1605d6ab1dc0SHong Zhang     }
1606d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1607d6ab1dc0SHong Zhang   }
1608d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1609d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1610d6ab1dc0SHong Zhang 
1611d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1612d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1613d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1614d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1615d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1616d6ab1dc0SHong Zhang }
1617d6ab1dc0SHong Zhang 
1618c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
16192bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
16202bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1621d6ab1dc0SHong Zhang #undef __FUNCT__
16226da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
16236da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1624d6ab1dc0SHong Zhang {
1625d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1626d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1627d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
16280298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1629d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1630d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1631d6ab1dc0SHong Zhang   PetscInt            nnz;
1632d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1633d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1634ce94432eSBarry Smith   MPI_Comm            comm;
1635d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1636d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1637d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1638d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1639d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1640d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1641d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1642d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1643d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1644d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
16454b38b95cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1646d6ab1dc0SHong Zhang   PetscScalar         *vals;
1647d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
16480acc65b4SHong Zhang   PetscTable          ta;
1649d6ab1dc0SHong Zhang 
1650d6ab1dc0SHong Zhang   PetscFunctionBegin;
1651ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1652d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
16536c4ed002SBarry Smith   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1654d6ab1dc0SHong Zhang 
1655d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1656d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1657d6ab1dc0SHong Zhang 
1658d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1659b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1660d6ab1dc0SHong Zhang 
1661d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1662d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16632205254eSKarl Rupp 
1664d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1665d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1666d6ab1dc0SHong Zhang   ai          = a_loc->i;
1667d6ab1dc0SHong Zhang   aj          = a_loc->j;
1668d6ab1dc0SHong Zhang 
1669d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1670d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1671d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1672d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1673d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1674d6ab1dc0SHong Zhang 
1675d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1676d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1677d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1678d6ab1dc0SHong Zhang 
1679d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1680d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1681d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1682854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1683d6ab1dc0SHong Zhang   coi[0] = 0;
1684d6ab1dc0SHong Zhang 
1685d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1686f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1687a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1688d6ab1dc0SHong Zhang   current_space = free_space;
168919f0e57aSHong Zhang 
1690d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1691c373ccc6SHong Zhang   ierr = PetscTableCreate(aN,aN,&ta);CHKERRQ(ierr);
16924b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
16930acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
16940acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1695d6ab1dc0SHong Zhang 
1696d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1697d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1698d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1699d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1700d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1701d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1702d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1703d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1704d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1705d6ab1dc0SHong Zhang     }
1706d6ab1dc0SHong Zhang     nnz = lnk[0];
1707d6ab1dc0SHong Zhang 
1708d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1709d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1710f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1711d6ab1dc0SHong Zhang       nspacedouble++;
1712d6ab1dc0SHong Zhang     }
1713d6ab1dc0SHong Zhang 
1714d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1715d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
17162205254eSKarl Rupp 
1717d6ab1dc0SHong Zhang     current_space->array           += nnz;
1718d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1719d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
17202205254eSKarl Rupp 
1721d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1722d6ab1dc0SHong Zhang   }
1723d6ab1dc0SHong Zhang 
1724854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1725d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
17260acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
17272205254eSKarl Rupp 
1728118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1729d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1730d6ab1dc0SHong Zhang 
1731d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1732d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1733d6ab1dc0SHong Zhang   /* determine row ownership */
1734b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1735d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17362205254eSKarl Rupp 
1737d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1738d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17392205254eSKarl Rupp 
1740d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1741d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1742d6ab1dc0SHong Zhang 
1743d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17441795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1745785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17462205254eSKarl Rupp 
1747d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1748d6ab1dc0SHong Zhang   merge->nsend = 0;
1749d6ab1dc0SHong Zhang 
1750854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1751d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1752d6ab1dc0SHong Zhang 
1753d6ab1dc0SHong Zhang   proc = 0;
1754d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1755d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1756d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1757d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1758d6ab1dc0SHong Zhang   }
1759d6ab1dc0SHong Zhang 
1760d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1761d6ab1dc0SHong Zhang   owners_co[0] = 0;
1762d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1763d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1764d6ab1dc0SHong Zhang     if (len_si[proc]) {
1765d6ab1dc0SHong Zhang       merge->nsend++;
1766d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1767d6ab1dc0SHong Zhang       len         += len_si[proc];
1768d6ab1dc0SHong Zhang     }
1769d6ab1dc0SHong Zhang   }
1770d6ab1dc0SHong Zhang 
1771d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17720298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1773d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1774d6ab1dc0SHong Zhang 
1775d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1776d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1777d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1778854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1779d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1780d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1781d6ab1dc0SHong Zhang     i    = owners_co[proc];
1782d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1783d6ab1dc0SHong Zhang     k++;
1784d6ab1dc0SHong Zhang   }
1785d6ab1dc0SHong Zhang 
1786d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1787785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1788d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1789c280ed6aSJed Brown     PetscMPIInt icompleted;
1790c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1791d6ab1dc0SHong Zhang   }
1792d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1793d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1794d6ab1dc0SHong Zhang 
17950acc65b4SHong Zhang   /* add received column indices into table to update Armax */
17960acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
17970acc65b4SHong Zhang     Jptr = buf_rj[k];
17980acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
17990acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
18000acc65b4SHong Zhang     }
18010acc65b4SHong Zhang   }
18020acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
18030acc65b4SHong Zhang 
1804d6ab1dc0SHong Zhang   /* send and recv coi */
1805d6ab1dc0SHong Zhang   /*-------------------*/
1806d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1807d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1808854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1809d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1810d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1811d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1812d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1813d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1814d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1815d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1816d6ab1dc0SHong Zhang     */
1817d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1818d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1819d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1820d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1821d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1822d6ab1dc0SHong Zhang     nrows       = 0;
1823d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1824d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1825d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1826d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1827d6ab1dc0SHong Zhang       nrows++;
1828d6ab1dc0SHong Zhang     }
1829d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1830d6ab1dc0SHong Zhang     k++;
1831d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1832d6ab1dc0SHong Zhang   }
1833d6ab1dc0SHong Zhang   i = merge->nrecv;
1834d6ab1dc0SHong Zhang   while (i--) {
1835c280ed6aSJed Brown     PetscMPIInt icompleted;
1836c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1837d6ab1dc0SHong Zhang   }
1838d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1839d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1840d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1841d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1842d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1843d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1844d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1845d6ab1dc0SHong Zhang 
1846d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1847d6ab1dc0SHong Zhang   /*------------------------------------------*/
1848d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1849854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1850d6ab1dc0SHong Zhang   bi[0] = 0;
1851d6ab1dc0SHong Zhang 
1852d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1853f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1854a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1855d6ab1dc0SHong Zhang   current_space = free_space;
1856d6ab1dc0SHong Zhang 
1857dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1858d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1859d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1860d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1861d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18622205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1863d6ab1dc0SHong Zhang   }
1864d6ab1dc0SHong Zhang 
18650acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1866d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1867d6ab1dc0SHong Zhang   rmax = 0;
1868d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1869d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1870d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1871d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1872d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1873d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1874d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1875d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1876d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1877d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1878d6ab1dc0SHong Zhang     }
1879d6ab1dc0SHong Zhang 
1880d6ab1dc0SHong Zhang     /* add received col data into lnk */
1881d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1882d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1883d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1884d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1885d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1886d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1887d6ab1dc0SHong Zhang       }
1888d6ab1dc0SHong Zhang     }
1889d6ab1dc0SHong Zhang     nnz = lnk[0];
1890d6ab1dc0SHong Zhang 
1891d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1892d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1893f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1894d6ab1dc0SHong Zhang       nspacedouble++;
1895d6ab1dc0SHong Zhang     }
1896d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1897d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1898d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18992205254eSKarl Rupp 
1900d6ab1dc0SHong Zhang     current_space->array           += nnz;
1901d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1902d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
19032205254eSKarl Rupp 
1904d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1905d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1906d6ab1dc0SHong Zhang   }
1907f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1908d6ab1dc0SHong Zhang 
1909854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1910d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1911118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1912d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1913d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
19140acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
19150acc65b4SHong Zhang 
1916d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1917d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1918d6ab1dc0SHong Zhang 
1919d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1920d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
19211795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1922d6ab1dc0SHong Zhang 
1923d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1924d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
192533d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1926d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1927d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1928d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1929d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1930d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1931d6ab1dc0SHong Zhang     row  = i + rstart;
1932d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1933d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1934d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1935d6ab1dc0SHong Zhang   }
1936d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1937d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1938d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1939d6ab1dc0SHong Zhang 
1940d6ab1dc0SHong Zhang   merge->bi        = bi;
1941d6ab1dc0SHong Zhang   merge->bj        = bj;
1942d6ab1dc0SHong Zhang   merge->coi       = coi;
1943d6ab1dc0SHong Zhang   merge->coj       = coj;
1944d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1945d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1946d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1947d6ab1dc0SHong Zhang 
1948d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1949d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19502205254eSKarl Rupp 
1951d6ab1dc0SHong Zhang   c->ptap     = ptap;
19520298fd71SBarry Smith   ptap->api   = NULL;
19530298fd71SBarry Smith   ptap->apj   = NULL;
1954d6ab1dc0SHong Zhang   ptap->merge = merge;
19550298fd71SBarry Smith   ptap->apa   = NULL;
1956a560ef98SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1957a560ef98SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1958a560ef98SHong Zhang 
1959a560ef98SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1960a560ef98SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1961a560ef98SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1962f96d37fcSHong Zhang 
1963f96d37fcSHong Zhang   *C = Cmpi;
1964f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1965f96d37fcSHong Zhang   if (bi[pn] != 0) {
196657622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
196757622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1968f96d37fcSHong Zhang   } else {
1969f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1970f96d37fcSHong Zhang   }
1971f96d37fcSHong Zhang #endif
1972f96d37fcSHong Zhang   PetscFunctionReturn(0);
1973f96d37fcSHong Zhang }
1974