xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision e55be12dc9957904bffc077bd8358204ed1851a4)
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 
132499ec78SHong Zhang #undef __FUNCT__
142499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
152499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
162499ec78SHong Zhang {
172499ec78SHong Zhang   PetscErrorCode ierr;
180fc8cf34SHong Zhang   const char     *algTypes[2] = {"scalable","nonscalable"};
190d3441aeSHong Zhang   PetscInt       alg=1; /* set default algorithm */
20*e55be12dSHong Zhang   MPI_Comm       comm;
212499ec78SHong Zhang 
222499ec78SHong Zhang   PetscFunctionBegin;
232499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
24*e55be12dSHong Zhang     ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
25*e55be12dSHong Zhang     if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) {
26*e55be12dSHong Zhang       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);
27*e55be12dSHong Zhang     }
28*e55be12dSHong Zhang 
290fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
3078d30b94SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr);
310fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
320fc8cf34SHong Zhang 
333ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
340fc8cf34SHong Zhang     switch (alg) {
350fc8cf34SHong Zhang     case 1:
360fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
370fc8cf34SHong Zhang       break;
380fc8cf34SHong Zhang     default:
39b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
400fc8cf34SHong Zhang       break;
410fc8cf34SHong Zhang     }
423ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
432499ec78SHong Zhang   }
443ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
45598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
463ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
472499ec78SHong Zhang   PetscFunctionReturn(0);
482499ec78SHong Zhang }
492499ec78SHong Zhang 
50f33d1a9aSHong Zhang #undef __FUNCT__
51a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
52a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
53598bc09dSHong Zhang {
54598bc09dSHong Zhang   PetscErrorCode ierr;
55598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
56598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
57598bc09dSHong Zhang 
58598bc09dSHong Zhang   PetscFunctionBegin;
59b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
60598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
61a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
62a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
63ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
64a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
65a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
66598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
67598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
68598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
69598bc09dSHong Zhang   PetscFunctionReturn(0);
70598bc09dSHong Zhang }
71598bc09dSHong Zhang 
722499ec78SHong Zhang #undef __FUNCT__
73a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
74a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
754ae0847bSHong Zhang {
764ae0847bSHong Zhang   PetscErrorCode ierr;
774ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
784ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
794ae0847bSHong Zhang 
804ae0847bSHong Zhang   PetscFunctionBegin;
814ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
8226fbe8dcSKarl Rupp 
834ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
844ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
854ae0847bSHong Zhang   PetscFunctionReturn(0);
864ae0847bSHong Zhang }
874ae0847bSHong Zhang 
884ae0847bSHong Zhang #undef __FUNCT__
89f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
90f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
91598bc09dSHong Zhang {
92598bc09dSHong Zhang   PetscErrorCode ierr;
934ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
94598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
95598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
96c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
97598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
989ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
99c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
100598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
101c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
10229825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
103598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1049467f45fSHong Zhang   MPI_Comm       comm;
1059467f45fSHong Zhang   PetscMPIInt    size;
106598bc09dSHong Zhang 
107598bc09dSHong Zhang   PetscFunctionBegin;
1089467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1099467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1109467f45fSHong Zhang 
111a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
112598bc09dSHong Zhang   /*-----------------------------------------------------*/
113a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
114b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
115a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
116598bc09dSHong Zhang 
117598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
118598bc09dSHong Zhang   /*----------------------------------------------------------*/
119598bc09dSHong Zhang   /* get data from symbolic products */
120a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
121c12557a1SHong Zhang   p_oth = NULL;
1229467f45fSHong Zhang   if (size >1) {
1239467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1249467f45fSHong Zhang   }
125598bc09dSHong Zhang 
126598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
127598bc09dSHong Zhang   apa = ptap->apa;
128598bc09dSHong Zhang 
12929825010SHong Zhang   api = ptap->api;
13029825010SHong Zhang   apj = ptap->apj;
131598bc09dSHong Zhang   for (i=0; i<cm; i++) {
132c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
133e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
134598bc09dSHong Zhang 
135598bc09dSHong Zhang     /* set values in C */
136598bc09dSHong Zhang     apJ  = apj + api[i];
137598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
138598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
139598bc09dSHong Zhang 
140598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
141598bc09dSHong Zhang     ca = coa + co->i[i];
142598bc09dSHong Zhang     k  = 0;
143598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
144598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
145598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
146598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
147598bc09dSHong Zhang       k++;
148598bc09dSHong Zhang     }
149598bc09dSHong Zhang 
150598bc09dSHong Zhang     /* diagonal part of C */
151598bc09dSHong Zhang     ca = cda + cd->i[i];
152598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
153598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
154598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
155598bc09dSHong Zhang       k++;
156598bc09dSHong Zhang     }
157598bc09dSHong Zhang 
158598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
159598bc09dSHong Zhang     ca = coa + co->i[i];
160598bc09dSHong Zhang     for (; k0<conz; k0++) {
161598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
162598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
163598bc09dSHong Zhang       k++;
164598bc09dSHong Zhang     }
165598bc09dSHong Zhang   }
166598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168598bc09dSHong Zhang   PetscFunctionReturn(0);
169598bc09dSHong Zhang }
170598bc09dSHong Zhang 
171598bc09dSHong Zhang #undef __FUNCT__
1720fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
1730fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
174598bc09dSHong Zhang {
175598bc09dSHong Zhang   PetscErrorCode     ierr;
176ce94432eSBarry Smith   MPI_Comm           comm;
1779467f45fSHong Zhang   PetscMPIInt        size;
178bfcd1a73SHong Zhang   Mat                Cmpi;
179598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1800298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1814ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
182bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1834ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
1844ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
185d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
186*e55be12dSHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n,Crmax;
187598bc09dSHong Zhang   PetscBT            lnkbt;
188598bc09dSHong Zhang   PetscScalar        *apa;
189bfcd1a73SHong Zhang   PetscReal          afill;
190*e55be12dSHong Zhang   PetscTable         ta;
191598bc09dSHong Zhang 
192598bc09dSHong Zhang   PetscFunctionBegin;
193ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1949467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1959467f45fSHong Zhang 
196598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
197b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
198598bc09dSHong Zhang 
199598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
200b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
20119f0e57aSHong Zhang 
202598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
203a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
204598bc09dSHong Zhang 
205a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
206598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2079467f45fSHong Zhang   if (size > 1) {
2089467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
209598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
21020e3dc75SHong Zhang   } else {
21120e3dc75SHong Zhang     p_oth = NULL;
21220e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2139467f45fSHong Zhang   }
214598bc09dSHong Zhang 
215598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
216598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
217854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
218a1a4e84aSHong Zhang   ptap->api = api;
219598bc09dSHong Zhang   api[0]    = 0;
220598bc09dSHong Zhang 
221598bc09dSHong Zhang   /* create and initialize a linked list */
222*e55be12dSHong Zhang   Crmax = p_loc->rmax+p_oth->rmax;
223*e55be12dSHong Zhang   ierr = PetscTableCreate(2*Crmax,pN,&ta);CHKERRQ(ierr);
224*e55be12dSHong Zhang   for (row=0; row<ptap->P_loc->rmap->N; row++) {
225*e55be12dSHong Zhang     nzi = p_loc->i[row+1] - p_loc->i[row];
226*e55be12dSHong Zhang     for (j=0; j<nzi; j++) {
227*e55be12dSHong Zhang       Jptr = j + p_loc->j + p_loc->i[row];
228*e55be12dSHong Zhang       ierr = PetscTableAdd(ta,*Jptr+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
2299467f45fSHong Zhang     }
230*e55be12dSHong Zhang   }
231*e55be12dSHong Zhang   for (row=0; row<ptap->P_oth->rmap->N; row++) {
232*e55be12dSHong Zhang     nzi = p_oth->i[row+1] - p_oth->i[row];
233*e55be12dSHong Zhang     for (j=0; j<nzi; j++) {
234*e55be12dSHong Zhang       Jptr = j + p_oth->j + p_oth->i[row];
235*e55be12dSHong Zhang       ierr = PetscTableAdd(ta,*Jptr+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
236*e55be12dSHong Zhang     }
237*e55be12dSHong Zhang   }
238*e55be12dSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
239*e55be12dSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
240*e55be12dSHong Zhang 
241*e55be12dSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
242598bc09dSHong Zhang 
243bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
244bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+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) {
272598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(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"
3479123193aSHong Zhang 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;
720f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
721cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
722cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
723cf1a0672SHong Zhang   PetscReal          afill;
724cf1a0672SHong Zhang   PetscScalar        *apa;
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 */
756cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
757a7c7454dSHong Zhang   if (size >1) {
758cf1a0672SHong Zhang     prmax = PetscMax(p_loc->rmax,p_oth->rmax);
759a7c7454dSHong Zhang   } else {
760a7c7454dSHong Zhang     prmax = p_loc->rmax;
761a7c7454dSHong Zhang   }
762cf1a0672SHong Zhang   nlnk_max = armax*prmax;
763cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
764f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
765cf1a0672SHong Zhang 
766cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
767cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
7682205254eSKarl Rupp 
769cf1a0672SHong Zhang   current_space = free_space;
770cf1a0672SHong Zhang 
771cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
772cf1a0672SHong Zhang   for (i=0; i<am; i++) {
773cf1a0672SHong Zhang     /* diagonal portion of A */
774cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
775cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
776cf1a0672SHong Zhang       row  = *adj++;
777cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
778cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
779cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
780f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
781cf1a0672SHong Zhang     }
782cf1a0672SHong Zhang     /* off-diagonal portion of A */
783cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
784cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
785cf1a0672SHong Zhang       row  = *aoj++;
786cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
787cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
788f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
789cf1a0672SHong Zhang     }
790cf1a0672SHong Zhang 
791f84c536bSHong Zhang     apnz     = *lnk;
792cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
793cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
794cf1a0672SHong Zhang 
795cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
796cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
797cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
798cf1a0672SHong Zhang       nspacedouble++;
799cf1a0672SHong Zhang     }
800cf1a0672SHong Zhang 
801cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
802f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
803cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
8042205254eSKarl Rupp 
805cf1a0672SHong Zhang     current_space->array           += apnz;
806cf1a0672SHong Zhang     current_space->local_used      += apnz;
807cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
808cf1a0672SHong Zhang   }
809cf1a0672SHong Zhang 
810cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
811cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
812854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
813cf1a0672SHong Zhang   apj  = ptap->apj;
814cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
815f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
816cf1a0672SHong Zhang 
817cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
818cf1a0672SHong Zhang   /*----------------------------------------------------*/
819cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
820cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
82133d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
822cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
823cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
824cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
825cf1a0672SHong Zhang 
82629825010SHong Zhang   /* malloc apa for assembly Cmpi */
8271795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8282205254eSKarl Rupp 
82929825010SHong Zhang   ptap->apa = apa;
830cf1a0672SHong Zhang   for (i=0; i<am; i++) {
831cf1a0672SHong Zhang     row  = i + rstart;
832cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
833cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
834cf1a0672SHong Zhang     apj += apnz;
835cf1a0672SHong Zhang   }
836cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
837cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
838cf1a0672SHong Zhang 
839cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
840cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
841f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
842cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
843cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
844cf1a0672SHong Zhang 
845cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
846cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
847cf1a0672SHong Zhang   c->ptap = ptap;
848cf1a0672SHong Zhang 
849cf1a0672SHong Zhang   *C = Cmpi;
850cf1a0672SHong Zhang 
851cf1a0672SHong Zhang   /* set MatInfo */
852118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
853cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
854cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
855cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
856cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
857cf1a0672SHong Zhang 
858cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
859cf1a0672SHong Zhang   if (api[am]) {
86057622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
86157622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
862cf1a0672SHong Zhang   } else {
863cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
864cf1a0672SHong Zhang   }
865cf1a0672SHong Zhang #endif
8661634b032SHong Zhang   PetscFunctionReturn(0);
8671634b032SHong Zhang }
868f96d37fcSHong Zhang 
869f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
870f96d37fcSHong Zhang #undef __FUNCT__
871f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
872c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
873f96d37fcSHong Zhang {
874f96d37fcSHong Zhang   PetscErrorCode ierr;
875c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
876c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
877f96d37fcSHong Zhang 
878f96d37fcSHong Zhang   PetscFunctionBegin;
879c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
880c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
881c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
882c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
883c216dbf3SHong Zhang 
884c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
885c216dbf3SHong Zhang     switch (alg) {
886c216dbf3SHong Zhang     case 1:
8872bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
888c216dbf3SHong Zhang       break;
889c216dbf3SHong Zhang     case 2:
890275476c6SMatthew G. Knepley     {
891ab1b013aSHong Zhang       Mat         Pt;
892ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
893ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
894ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
895ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
896ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
897ab1b013aSHong Zhang       ptap     = c->ptap;
898ab1b013aSHong Zhang       ptap->Pt = Pt;
899c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
900c216dbf3SHong Zhang       PetscFunctionReturn(0);
901275476c6SMatthew G. Knepley     }
902c216dbf3SHong Zhang       break;
903c216dbf3SHong Zhang     default:
9046da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
905c216dbf3SHong Zhang       break;
906c216dbf3SHong Zhang     }
907c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
908c216dbf3SHong Zhang   }
909c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
910c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
911c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
912ab1b013aSHong Zhang   PetscFunctionReturn(0);
913ab1b013aSHong Zhang }
914ab1b013aSHong Zhang 
915c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
916c216dbf3SHong Zhang #undef __FUNCT__
917c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
918c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
919c216dbf3SHong Zhang {
920c216dbf3SHong Zhang   PetscErrorCode ierr;
9212bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9222bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9232bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
924c216dbf3SHong Zhang 
925c216dbf3SHong Zhang   PetscFunctionBegin;
926c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
927c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
928f96d37fcSHong Zhang   PetscFunctionReturn(0);
929f96d37fcSHong Zhang }
930f96d37fcSHong Zhang 
9316da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
932f96d37fcSHong Zhang #undef __FUNCT__
9336da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9346da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
935f96d37fcSHong Zhang {
936c5af039cSHong Zhang   PetscErrorCode      ierr;
937c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
938497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
939c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
940c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
941d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
942497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
943e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
944c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
945ce94432eSBarry Smith   MPI_Comm            comm;
946c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
947c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
948c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
949c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
950c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
951c5af039cSHong Zhang   MPI_Status          *status;
952c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
953d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
954a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
955c5af039cSHong Zhang   Mat                 A_loc;
956c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
957f96d37fcSHong Zhang 
958f96d37fcSHong Zhang   PetscFunctionBegin;
959ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
960c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
961c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
962c5af039cSHong Zhang 
963c5af039cSHong Zhang   ptap  = c->ptap;
964c5af039cSHong Zhang   merge = ptap->merge;
965c5af039cSHong Zhang 
966c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
967c5af039cSHong Zhang   /*--------------------------------------------------------------*/
968c5af039cSHong Zhang   /* get data from symbolic products */
969c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
970854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
971c5af039cSHong Zhang 
972c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
973c5af039cSHong Zhang   owners = merge->rowmap->range;
974854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
975c5af039cSHong Zhang 
976c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
977c5af039cSHong Zhang   A_loc = ptap->A_loc;
978c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
979c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
980d6ab1dc0SHong Zhang   ai    = a_loc->i;
981d6ab1dc0SHong Zhang   aj    = a_loc->j;
982c5af039cSHong Zhang 
983854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
984c5af039cSHong Zhang 
985c5af039cSHong Zhang   for (i=0; i<am; i++) {
986e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
987d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
988d6ab1dc0SHong Zhang     adj = aj + ai[i];
989d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
990c5af039cSHong Zhang     for (j=0; j<anz; j++) {
991e745005bSHong Zhang       aval[adj[j]] = ada[j];
992c5af039cSHong Zhang     }
993c5af039cSHong Zhang 
994c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
995c5af039cSHong Zhang     /*--------------------------------------------------------------*/
996c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
997c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
998c5af039cSHong Zhang     poJ = po->j + po->i[i];
999c5af039cSHong Zhang     pA  = po->a + po->i[i];
1000c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1001c5af039cSHong Zhang       row = poJ[j];
1002c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
1003c5af039cSHong Zhang       cj  = coj + coi[row];
1004c5af039cSHong Zhang       ca  = coa + coi[row];
1005c5af039cSHong Zhang       /* perform dense axpy */
1006c5af039cSHong Zhang       valtmp = pA[j];
1007c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1008e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1009c5af039cSHong Zhang       }
1010c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1011c5af039cSHong Zhang     }
1012c5af039cSHong Zhang 
1013c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
1014c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1015c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1016c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1017c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1018c5af039cSHong Zhang       row = pdJ[j];
1019c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1020c5af039cSHong Zhang       cj  = bj + bi[row];
1021c5af039cSHong Zhang       ca  = ba + bi[row];
1022c5af039cSHong Zhang       /* perform dense axpy */
1023c5af039cSHong Zhang       valtmp = pA[j];
1024c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1025e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1026c5af039cSHong Zhang       }
1027c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1028c5af039cSHong Zhang     }
1029c5af039cSHong Zhang 
1030d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1031d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1032e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1033c5af039cSHong Zhang   }
1034c5af039cSHong Zhang 
1035c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1036c5af039cSHong Zhang   /*------------------------------------*/
1037c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1038c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1039c5af039cSHong Zhang   len_s  = merge->len_s;
1040c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1041c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1042c5af039cSHong Zhang 
1043dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1044c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1045c5af039cSHong Zhang     if (!len_s[proc]) continue;
1046c5af039cSHong Zhang     i    = merge->owners_co[proc];
1047c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1048c5af039cSHong Zhang     k++;
1049c5af039cSHong Zhang   }
1050c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1051c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1052c5af039cSHong Zhang 
1053c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1054c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1055c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1056c5af039cSHong Zhang 
1057c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1058c5af039cSHong Zhang   /*----------------------------------------------------*/
1059dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1060c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1061c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1062c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1063c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1064c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1065c5af039cSHong Zhang   }
1066c5af039cSHong Zhang 
1067c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1068c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1069c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1070c5af039cSHong Zhang     ba_i = ba + bi[i];
1071c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1072c5af039cSHong Zhang     /* add received vals into ba */
1073c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1074c5af039cSHong Zhang       /* i-th row */
1075c5af039cSHong Zhang       if (i == *nextrow[k]) {
1076c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1077c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1078c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1079c5af039cSHong Zhang         nextcj = 0;
1080c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1081c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1082c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1083c5af039cSHong Zhang           }
1084c5af039cSHong Zhang         }
1085c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1086c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1087c5af039cSHong Zhang       }
1088c5af039cSHong Zhang     }
1089c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1090c5af039cSHong Zhang   }
1091c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1092c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1093c5af039cSHong Zhang 
1094c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1095c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1096c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1097c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1098e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1099f96d37fcSHong Zhang   PetscFunctionReturn(0);
1100f96d37fcSHong Zhang }
1101f96d37fcSHong Zhang 
1102c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1103c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1104f96d37fcSHong Zhang #undef __FUNCT__
11052bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
11062bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1107f96d37fcSHong Zhang {
1108f96d37fcSHong Zhang   PetscErrorCode      ierr;
11094e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1110f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
11110298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1112f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1113f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1114f96d37fcSHong Zhang   PetscInt            nnz;
11154e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1116497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1117f96d37fcSHong Zhang   PetscBT             lnkbt;
1118ce94432eSBarry Smith   MPI_Comm            comm;
1119f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1120f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1121f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1122f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1123f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1124f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1125f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1126f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1127d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1128f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1129438d860cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1130f96d37fcSHong Zhang   PetscScalar         *vals;
11314e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1132f96d37fcSHong Zhang 
1133f96d37fcSHong Zhang   PetscFunctionBegin;
1134ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1135c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1136c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1137c5af039cSHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1138c5af039cSHong Zhang   }
1139c5af039cSHong Zhang 
1140f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1141f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1142f96d37fcSHong Zhang 
1143f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1144b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1145f96d37fcSHong Zhang 
1146f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1147f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11482205254eSKarl Rupp 
1149c5af039cSHong Zhang   ptap->A_loc = A_loc;
11502205254eSKarl Rupp 
11511c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1152d6ab1dc0SHong Zhang   ai    = a_loc->i;
1153d6ab1dc0SHong Zhang   aj    = a_loc->j;
1154f96d37fcSHong Zhang 
1155f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1156f96d37fcSHong Zhang   /*----------------------------------------------------*/
11574e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11584e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11594e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11604e938277SHong Zhang 
11614e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11624e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11634e938277SHong Zhang   poti = pot->i; potj = pot->j;
1164f96d37fcSHong Zhang 
1165f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11662205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1167854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1168f96d37fcSHong Zhang   coi[0] = 0;
1169f96d37fcSHong Zhang 
1170f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1171d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1172a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1173f96d37fcSHong Zhang   current_space = free_space;
117419f0e57aSHong Zhang 
1175f96d37fcSHong Zhang   /* create and initialize a linked list */
1176438d860cSHong Zhang   ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1177f96d37fcSHong Zhang 
1178f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1179f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1180f96d37fcSHong Zhang     ptJ = potj + poti[i];
1181f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1182f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1183d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1184d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1185f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1186d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1187f96d37fcSHong Zhang     }
11884e938277SHong Zhang     nnz = lnk[0];
1189f96d37fcSHong Zhang 
1190f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1191f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1192f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1193f96d37fcSHong Zhang       nspacedouble++;
1194f96d37fcSHong Zhang     }
1195f96d37fcSHong Zhang 
1196f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11974e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11982205254eSKarl Rupp 
1199f96d37fcSHong Zhang     current_space->array           += nnz;
1200f96d37fcSHong Zhang     current_space->local_used      += nnz;
1201f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
12022205254eSKarl Rupp 
1203f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1204f96d37fcSHong Zhang   }
1205f96d37fcSHong Zhang 
1206854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1207f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
12082205254eSKarl Rupp 
1209118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1210f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1211f96d37fcSHong Zhang 
1212f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1213f96d37fcSHong Zhang   /*----------------------------------------------*/
1214f96d37fcSHong Zhang   /* determine row ownership */
1215b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1216f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12172205254eSKarl Rupp 
1218f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1219f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12202205254eSKarl Rupp 
1221f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1222f96d37fcSHong Zhang   owners = merge->rowmap->range;
1223f96d37fcSHong Zhang 
1224f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12251795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1226785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12272205254eSKarl Rupp 
1228f96d37fcSHong Zhang   len_s        = merge->len_s;
1229f96d37fcSHong Zhang   merge->nsend = 0;
1230f96d37fcSHong Zhang 
1231854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1232f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1233f96d37fcSHong Zhang 
1234f96d37fcSHong Zhang   proc = 0;
1235f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1236f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1237f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1238f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1239f96d37fcSHong Zhang   }
1240f96d37fcSHong Zhang 
1241f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1242f96d37fcSHong Zhang   owners_co[0] = 0;
1243f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1244f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1245f96d37fcSHong Zhang     if (len_si[proc]) {
1246f96d37fcSHong Zhang       merge->nsend++;
1247f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1248f96d37fcSHong Zhang       len         += len_si[proc];
1249f96d37fcSHong Zhang     }
1250f96d37fcSHong Zhang   }
1251f96d37fcSHong Zhang 
1252f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12530298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1254f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1255f96d37fcSHong Zhang 
1256f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1257f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1258f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1259854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1260f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1261f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1262f96d37fcSHong Zhang     i    = owners_co[proc];
1263f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1264f96d37fcSHong Zhang     k++;
1265f96d37fcSHong Zhang   }
1266f96d37fcSHong Zhang 
1267f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1268785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1269f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1270c280ed6aSJed Brown     PetscMPIInt icompleted;
1271c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1272f96d37fcSHong Zhang   }
1273f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1274f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1275f96d37fcSHong Zhang 
1276f96d37fcSHong Zhang   /* send and recv coi */
1277f96d37fcSHong Zhang   /*-------------------*/
1278f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1279f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1280854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1281f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1282f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1283f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1284f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1285f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1286f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1287f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1288f96d37fcSHong Zhang     */
1289f96d37fcSHong Zhang     /*-------------------------------------------*/
1290f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1291f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1292f96d37fcSHong Zhang     buf_si[0]   = nrows;
1293f96d37fcSHong Zhang     buf_si_i[0] = 0;
1294f96d37fcSHong Zhang     nrows       = 0;
1295f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1296f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1297f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1298f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1299f96d37fcSHong Zhang       nrows++;
1300f96d37fcSHong Zhang     }
1301f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1302f96d37fcSHong Zhang     k++;
1303f96d37fcSHong Zhang     buf_si += len_si[proc];
1304f96d37fcSHong Zhang   }
1305f96d37fcSHong Zhang   i = merge->nrecv;
1306f96d37fcSHong Zhang   while (i--) {
1307c280ed6aSJed Brown     PetscMPIInt icompleted;
1308c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1309f96d37fcSHong Zhang   }
1310f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1311f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1312f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1313f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1314f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1315f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1316f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1317f96d37fcSHong Zhang 
1318f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1319f96d37fcSHong Zhang   /*------------------------------------------*/
1320f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1321854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1322f96d37fcSHong Zhang   bi[0] = 0;
1323f96d37fcSHong Zhang 
1324c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1325d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1326a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1327f96d37fcSHong Zhang   current_space = free_space;
1328f96d37fcSHong Zhang 
1329dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1330f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1331f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1332f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1333f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1334f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1335f96d37fcSHong Zhang   }
1336f96d37fcSHong Zhang 
13371c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1338f96d37fcSHong Zhang   rmax = 0;
1339f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1340f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1341f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1342f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1343f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1344f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1345d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1346d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1347f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1348d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1349f96d37fcSHong Zhang     }
13504e938277SHong Zhang 
1351f96d37fcSHong Zhang     /* add received col data into lnk */
1352f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1353f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1354f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1355f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13564e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1357f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1358f96d37fcSHong Zhang       }
1359f96d37fcSHong Zhang     }
13604e938277SHong Zhang     nnz = lnk[0];
1361f96d37fcSHong Zhang 
1362f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1363f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1364f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1365f96d37fcSHong Zhang       nspacedouble++;
1366f96d37fcSHong Zhang     }
1367f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13684e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1369f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13702205254eSKarl Rupp 
1371f96d37fcSHong Zhang     current_space->array           += nnz;
1372f96d37fcSHong Zhang     current_space->local_used      += nnz;
1373f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13742205254eSKarl Rupp 
1375f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1376f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1377f96d37fcSHong Zhang   }
1378f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1379f96d37fcSHong Zhang 
1380854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1381f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13822205254eSKarl Rupp 
1383118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1384f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1385d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13864e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13874e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1388f96d37fcSHong Zhang 
13891c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13901c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13911795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1392f96d37fcSHong Zhang 
1393f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1394f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
139533d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1396f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1397f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1398f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1399f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1400f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1401f96d37fcSHong Zhang     row  = i + rstart;
1402f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1403f96d37fcSHong Zhang     Jptr = bj + bi[i];
1404f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1405f96d37fcSHong Zhang   }
1406f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1407f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1408f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1409f96d37fcSHong Zhang 
1410f96d37fcSHong Zhang   merge->bi        = bi;
1411f96d37fcSHong Zhang   merge->bj        = bj;
1412f96d37fcSHong Zhang   merge->coi       = coi;
1413f96d37fcSHong Zhang   merge->coj       = coj;
1414f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1415f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1416f96d37fcSHong Zhang   merge->owners_co = owners_co;
1417f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1418f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1419f96d37fcSHong Zhang 
14206da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1421c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1422c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1423f96d37fcSHong Zhang 
1424f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1425f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1426f96d37fcSHong Zhang   c->ptap     = ptap;
14270298fd71SBarry Smith   ptap->api   = NULL;
14280298fd71SBarry Smith   ptap->apj   = NULL;
1429f96d37fcSHong Zhang   ptap->merge = merge;
1430d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1431d6ab1dc0SHong Zhang 
1432d6ab1dc0SHong Zhang   *C = Cmpi;
1433d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1434d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
143557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
143657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1437d6ab1dc0SHong Zhang   } else {
1438d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1439d6ab1dc0SHong Zhang   }
1440d6ab1dc0SHong Zhang #endif
1441d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1442d6ab1dc0SHong Zhang }
1443d6ab1dc0SHong Zhang 
1444d6ab1dc0SHong Zhang #undef __FUNCT__
14456da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14466da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1447d6ab1dc0SHong Zhang {
1448d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1449d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1450d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1451d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1452d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1453e745005bSHong Zhang   PetscInt            *adj;
1454e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1455e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1456d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1457ce94432eSBarry Smith   MPI_Comm            comm;
1458d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1459d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1460d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1461d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1462d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1463d6ab1dc0SHong Zhang   MPI_Status          *status;
1464d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1465d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1466a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1467d6ab1dc0SHong Zhang   Mat                 A_loc;
1468d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1469d6ab1dc0SHong Zhang 
1470d6ab1dc0SHong Zhang   PetscFunctionBegin;
1471ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1472d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1473d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1474d6ab1dc0SHong Zhang 
1475d6ab1dc0SHong Zhang   ptap  = c->ptap;
1476d6ab1dc0SHong Zhang   merge = ptap->merge;
1477d6ab1dc0SHong Zhang 
1478e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1479e745005bSHong Zhang   /*------------------------------------------*/
1480d6ab1dc0SHong Zhang   /* get data from symbolic products */
1481d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1482854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1483d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1484d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14851795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1486d6ab1dc0SHong Zhang 
1487d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1488d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1489d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1490d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1491d6ab1dc0SHong Zhang   ai    = a_loc->i;
1492d6ab1dc0SHong Zhang   aj    = a_loc->j;
1493d6ab1dc0SHong Zhang 
1494d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1495d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1496d6ab1dc0SHong Zhang     adj = aj + ai[i];
1497d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1498d6ab1dc0SHong Zhang 
1499d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1500e745005bSHong Zhang     /*-------------------------------------------------------------*/
1501d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1502d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1503d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1504d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1505d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1506d6ab1dc0SHong Zhang       row = poJ[j];
1507d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1508d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1509e745005bSHong Zhang       /* perform sparse axpy */
1510e745005bSHong Zhang       nexta  = 0;
1511d6ab1dc0SHong Zhang       valtmp = pA[j];
1512e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1513e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1514e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1515e745005bSHong Zhang           nexta++;
1516d6ab1dc0SHong Zhang         }
1517e745005bSHong Zhang       }
1518e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1519d6ab1dc0SHong Zhang     }
1520d6ab1dc0SHong Zhang 
1521d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1522d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1523d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1524d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1525d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1526d6ab1dc0SHong Zhang       row = pdJ[j];
1527d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1528d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1529e745005bSHong Zhang       /* perform sparse axpy */
1530e745005bSHong Zhang       nexta  = 0;
1531d6ab1dc0SHong Zhang       valtmp = pA[j];
1532e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1533e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1534e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1535e745005bSHong Zhang           nexta++;
1536d6ab1dc0SHong Zhang         }
1537e745005bSHong Zhang       }
1538e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1539d6ab1dc0SHong Zhang     }
1540d6ab1dc0SHong Zhang   }
1541d6ab1dc0SHong Zhang 
1542d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1543d6ab1dc0SHong Zhang   /*------------------------------------*/
1544d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1545d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1546d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1547d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1548d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1549d6ab1dc0SHong Zhang 
1550dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1551d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1552d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1553d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1554d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1555d6ab1dc0SHong Zhang     k++;
1556d6ab1dc0SHong Zhang   }
1557d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1558d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1559d6ab1dc0SHong Zhang 
1560d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1561d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1562d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1563d6ab1dc0SHong Zhang 
1564d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1565d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1566dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1567d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1568e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1569d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1570d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1571d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1572d6ab1dc0SHong Zhang   }
1573d6ab1dc0SHong Zhang 
1574d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1575d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1576d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1577d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1578d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1579d6ab1dc0SHong Zhang     /* add received vals into ba */
1580d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1581d6ab1dc0SHong Zhang       /* i-th row */
1582d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1583d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1584d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1585d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1586d6ab1dc0SHong Zhang         nextcj = 0;
1587d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1588d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1589d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1590d6ab1dc0SHong Zhang           }
1591d6ab1dc0SHong Zhang         }
1592d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1593d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1594d6ab1dc0SHong Zhang       }
1595d6ab1dc0SHong Zhang     }
1596d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1597d6ab1dc0SHong Zhang   }
1598d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1599d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1600d6ab1dc0SHong Zhang 
1601d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1602d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1603d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1604d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1605d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1606d6ab1dc0SHong Zhang }
1607d6ab1dc0SHong Zhang 
16080acc65b4SHong Zhang #if 0
16090acc65b4SHong Zhang #undef __FUNCT__
16100acc65b4SHong Zhang #define __FUNCT__ "MatRowMergeMax_MPIAIJ"
16110acc65b4SHong Zhang PetscErrorCode MatRowMergeMax_MPIAIJ(Mat A,const PetscInt nexpect,PetscInt maxkey,PetscInt *rmax)
16120acc65b4SHong Zhang {
16130acc65b4SHong Zhang   PetscErrorCode  ierr;
16140acc65b4SHong Zhang   Mat_MPIAIJ      *a = (Mat_MPIAIJ*)A->data;
16150acc65b4SHong Zhang   Mat_SeqAIJ      *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
16160acc65b4SHong Zhang   PetscInt        am=A->rmap->n,*col,row,j,nz,*garray=a->garray;
16170acc65b4SHong Zhang   PetscTable      ta;
16180acc65b4SHong Zhang 
16190acc65b4SHong Zhang   PetscFunctionBegin;
16200acc65b4SHong Zhang   ierr = PetscTableCreate(nexpect,maxkey,&ta);CHKERRQ(ierr);
16210acc65b4SHong Zhang   /* diagonal part */
16220acc65b4SHong Zhang   for (row=0; row<am; row++) {
16230acc65b4SHong Zhang     nz = ad->i[row+1] - ad->i[row];
16240acc65b4SHong Zhang     for (j=0; j<nz; j++) {
16250acc65b4SHong Zhang       col = j + ad->j + ad->i[row];
16260acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*col+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
16270acc65b4SHong Zhang     }
16280acc65b4SHong Zhang   }
16290acc65b4SHong Zhang 
16300acc65b4SHong Zhang   /* off-diagonal part */
16310acc65b4SHong Zhang   for (row=0; row<am; row++) {
16320acc65b4SHong Zhang     nz = ao->i[row+1] - ao->i[row];
16330acc65b4SHong Zhang     for (j=0; j<nz; j++) {
16340acc65b4SHong Zhang       col = j + ao->j + ao->i[row];
16350acc65b4SHong Zhang       ierr = PetscTableAdd(ta,garray[*col]+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
16360acc65b4SHong Zhang     }
16370acc65b4SHong Zhang   }
16380acc65b4SHong Zhang 
16390acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,rmax);CHKERRQ(ierr);
16400acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
16410acc65b4SHong Zhang   printf("MatRowMergeMax_MPIAIJ...rmax %d\n",*rmax);
16420acc65b4SHong Zhang   PetscFunctionReturn(0);
16430acc65b4SHong Zhang }
16440acc65b4SHong Zhang #endif
16450acc65b4SHong Zhang 
1646c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
16472bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
16482bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1649d6ab1dc0SHong Zhang #undef __FUNCT__
16506da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
16516da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1652d6ab1dc0SHong Zhang {
1653d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1654d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1655d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
16560298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1657d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1658d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1659d6ab1dc0SHong Zhang   PetscInt            nnz;
1660d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1661d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1662ce94432eSBarry Smith   MPI_Comm            comm;
1663d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1664d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1665d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1666d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1667d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1668d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1669d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1670d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1671d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1672d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
16730acc65b4SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax,*col,nz;
1674d6ab1dc0SHong Zhang   PetscScalar         *vals;
1675d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
16760acc65b4SHong Zhang   PetscTable          ta;
1677d6ab1dc0SHong Zhang 
1678d6ab1dc0SHong Zhang   PetscFunctionBegin;
1679ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1680d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1681d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1682d6ab1dc0SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1683d6ab1dc0SHong Zhang   }
1684d6ab1dc0SHong Zhang 
1685d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1686d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1687d6ab1dc0SHong Zhang 
1688d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1689b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1690d6ab1dc0SHong Zhang 
1691d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1692d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16932205254eSKarl Rupp 
1694d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1695d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1696d6ab1dc0SHong Zhang   ai          = a_loc->i;
1697d6ab1dc0SHong Zhang   aj          = a_loc->j;
1698d6ab1dc0SHong Zhang 
1699d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1700d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1701d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1702d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1703d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1704d6ab1dc0SHong Zhang 
1705d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1706d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1707d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1708d6ab1dc0SHong Zhang 
1709d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1710d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1711d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1712854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1713d6ab1dc0SHong Zhang   coi[0] = 0;
1714d6ab1dc0SHong Zhang 
1715d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1716d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1717a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1718d6ab1dc0SHong Zhang   current_space = free_space;
171919f0e57aSHong Zhang 
1720d6ab1dc0SHong Zhang   /* create and initialize a linked list */
17210acc65b4SHong Zhang   ierr = PetscTableCreate(2*a_loc->rmax,aN,&ta);CHKERRQ(ierr);
17220acc65b4SHong Zhang   for (row=0; row<am; row++) {
17230acc65b4SHong Zhang     nz = ai[row+1] - ai[row];
17240acc65b4SHong Zhang     for (j=0; j<nz; j++) {
17250acc65b4SHong Zhang       col = j + aj + ai[row];
17260acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*col+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
17270acc65b4SHong Zhang     }
17280acc65b4SHong Zhang   }
17290acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
17300acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1731d6ab1dc0SHong Zhang 
1732d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1733d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1734d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1735d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1736d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1737d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1738d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1739d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1740d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1741d6ab1dc0SHong Zhang     }
1742d6ab1dc0SHong Zhang     nnz = lnk[0];
1743d6ab1dc0SHong Zhang 
1744d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1745d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1746d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1747d6ab1dc0SHong Zhang       nspacedouble++;
1748d6ab1dc0SHong Zhang     }
1749d6ab1dc0SHong Zhang 
1750d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1751d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
17522205254eSKarl Rupp 
1753d6ab1dc0SHong Zhang     current_space->array           += nnz;
1754d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1755d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
17562205254eSKarl Rupp 
1757d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1758d6ab1dc0SHong Zhang   }
1759d6ab1dc0SHong Zhang 
1760854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1761d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
17620acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
17632205254eSKarl Rupp 
1764118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1765d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1766d6ab1dc0SHong Zhang 
1767d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1768d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1769d6ab1dc0SHong Zhang   /* determine row ownership */
1770b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1771d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17722205254eSKarl Rupp 
1773d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1774d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17752205254eSKarl Rupp 
1776d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1777d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1778d6ab1dc0SHong Zhang 
1779d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17801795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1781785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17822205254eSKarl Rupp 
1783d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1784d6ab1dc0SHong Zhang   merge->nsend = 0;
1785d6ab1dc0SHong Zhang 
1786854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1787d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1788d6ab1dc0SHong Zhang 
1789d6ab1dc0SHong Zhang   proc = 0;
1790d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1791d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1792d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1793d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1794d6ab1dc0SHong Zhang   }
1795d6ab1dc0SHong Zhang 
1796d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1797d6ab1dc0SHong Zhang   owners_co[0] = 0;
1798d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1799d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1800d6ab1dc0SHong Zhang     if (len_si[proc]) {
1801d6ab1dc0SHong Zhang       merge->nsend++;
1802d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1803d6ab1dc0SHong Zhang       len         += len_si[proc];
1804d6ab1dc0SHong Zhang     }
1805d6ab1dc0SHong Zhang   }
1806d6ab1dc0SHong Zhang 
1807d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
18080298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1809d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1810d6ab1dc0SHong Zhang 
1811d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1812d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1813d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1814854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1815d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1816d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1817d6ab1dc0SHong Zhang     i    = owners_co[proc];
1818d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang     k++;
1820d6ab1dc0SHong Zhang   }
1821d6ab1dc0SHong Zhang 
1822d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1823785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1824d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1825c280ed6aSJed Brown     PetscMPIInt icompleted;
1826c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1827d6ab1dc0SHong Zhang   }
1828d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1829d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1830d6ab1dc0SHong Zhang 
18310acc65b4SHong Zhang   /* add received column indices into table to update Armax */
18320acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
18330acc65b4SHong Zhang     Jptr = buf_rj[k];
18340acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
18350acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
18360acc65b4SHong Zhang     }
18370acc65b4SHong Zhang   }
18380acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
18390acc65b4SHong Zhang 
1840d6ab1dc0SHong Zhang   /* send and recv coi */
1841d6ab1dc0SHong Zhang   /*-------------------*/
1842d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1843d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1844854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1845d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1846d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1847d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1848d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1849d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1850d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1851d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1852d6ab1dc0SHong Zhang     */
1853d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1854d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1855d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1856d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1857d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1858d6ab1dc0SHong Zhang     nrows       = 0;
1859d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1860d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1861d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1862d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1863d6ab1dc0SHong Zhang       nrows++;
1864d6ab1dc0SHong Zhang     }
1865d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1866d6ab1dc0SHong Zhang     k++;
1867d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1868d6ab1dc0SHong Zhang   }
1869d6ab1dc0SHong Zhang   i = merge->nrecv;
1870d6ab1dc0SHong Zhang   while (i--) {
1871c280ed6aSJed Brown     PetscMPIInt icompleted;
1872c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1873d6ab1dc0SHong Zhang   }
1874d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1875d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1876d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1877d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1878d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1879d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1880d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1881d6ab1dc0SHong Zhang 
1882d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1883d6ab1dc0SHong Zhang   /*------------------------------------------*/
1884d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1885854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1886d6ab1dc0SHong Zhang   bi[0] = 0;
1887d6ab1dc0SHong Zhang 
1888d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1889d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1890a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1891d6ab1dc0SHong Zhang   current_space = free_space;
1892d6ab1dc0SHong Zhang 
1893dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1894d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1895d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1896d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1897d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18982205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1899d6ab1dc0SHong Zhang   }
1900d6ab1dc0SHong Zhang 
19010acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1902d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1903d6ab1dc0SHong Zhang   rmax = 0;
1904d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1905d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1906d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1907d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1908d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1909d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1910d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1911d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1912d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1913d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1914d6ab1dc0SHong Zhang     }
1915d6ab1dc0SHong Zhang 
1916d6ab1dc0SHong Zhang     /* add received col data into lnk */
1917d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1918d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1919d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1920d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1921d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1922d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1923d6ab1dc0SHong Zhang       }
1924d6ab1dc0SHong Zhang     }
1925d6ab1dc0SHong Zhang     nnz = lnk[0];
1926d6ab1dc0SHong Zhang 
1927d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1928d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1929d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1930d6ab1dc0SHong Zhang       nspacedouble++;
1931d6ab1dc0SHong Zhang     }
1932d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1933d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1934d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
19352205254eSKarl Rupp 
1936d6ab1dc0SHong Zhang     current_space->array           += nnz;
1937d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1938d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
19392205254eSKarl Rupp 
1940d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1941d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1942d6ab1dc0SHong Zhang   }
1943f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1944d6ab1dc0SHong Zhang 
1945854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1946d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1947118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1948d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1949d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
19500acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
19510acc65b4SHong Zhang 
1952d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1953d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1954d6ab1dc0SHong Zhang 
1955d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1956d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
19571795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1958d6ab1dc0SHong Zhang 
1959d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1960d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
196133d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1962d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1963d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1964d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1965d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1966d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1967d6ab1dc0SHong Zhang     row  = i + rstart;
1968d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1969d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1970d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1971d6ab1dc0SHong Zhang   }
1972d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1973d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1974d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1975d6ab1dc0SHong Zhang 
1976d6ab1dc0SHong Zhang   merge->bi        = bi;
1977d6ab1dc0SHong Zhang   merge->bj        = bj;
1978d6ab1dc0SHong Zhang   merge->coi       = coi;
1979d6ab1dc0SHong Zhang   merge->coj       = coj;
1980d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1981d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1982d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1983d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1984d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1985d6ab1dc0SHong Zhang 
19866da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1987d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1988c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1989d6ab1dc0SHong Zhang 
1990d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1991d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19922205254eSKarl Rupp 
1993d6ab1dc0SHong Zhang   c->ptap     = ptap;
19940298fd71SBarry Smith   ptap->api   = NULL;
19950298fd71SBarry Smith   ptap->apj   = NULL;
1996d6ab1dc0SHong Zhang   ptap->merge = merge;
1997d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
19980298fd71SBarry Smith   ptap->apa   = NULL;
1999f96d37fcSHong Zhang 
2000f96d37fcSHong Zhang   *C = Cmpi;
2001f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
2002f96d37fcSHong Zhang   if (bi[pn] != 0) {
200357622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
200457622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
2005f96d37fcSHong Zhang   } else {
2006f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
2007f96d37fcSHong Zhang   }
2008f96d37fcSHong Zhang #endif
2009f96d37fcSHong Zhang   PetscFunctionReturn(0);
2010f96d37fcSHong Zhang }
2011