xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 9467f45fc2f06212a563b445f167c8c6a3337d2e)
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>
119a6d0b0bSJed Brown #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"};
190fc8cf34SHong Zhang   PetscInt       alg=0; /* set default algorithm */
202499ec78SHong Zhang 
212499ec78SHong Zhang   PetscFunctionBegin;
222499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
230fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
240fc8cf34SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr);
250fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
260fc8cf34SHong Zhang 
273ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
280fc8cf34SHong Zhang     switch (alg) {
290fc8cf34SHong Zhang     case 1:
300fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
310fc8cf34SHong Zhang       break;
320fc8cf34SHong Zhang     default:
33b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
340fc8cf34SHong Zhang       break;
350fc8cf34SHong Zhang     }
363ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
372499ec78SHong Zhang   }
383ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
39598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
403ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
412499ec78SHong Zhang   PetscFunctionReturn(0);
422499ec78SHong Zhang }
432499ec78SHong Zhang 
44f33d1a9aSHong Zhang #undef __FUNCT__
45a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
46a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
47598bc09dSHong Zhang {
48598bc09dSHong Zhang   PetscErrorCode ierr;
49598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
50598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
51598bc09dSHong Zhang 
52598bc09dSHong Zhang   PetscFunctionBegin;
53b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
54598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
55a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
56a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
57ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
58a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
59a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
60598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
61598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
62598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
63598bc09dSHong Zhang   PetscFunctionReturn(0);
64598bc09dSHong Zhang }
65598bc09dSHong Zhang 
662499ec78SHong Zhang #undef __FUNCT__
67a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
68a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
694ae0847bSHong Zhang {
704ae0847bSHong Zhang   PetscErrorCode ierr;
714ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
724ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
734ae0847bSHong Zhang 
744ae0847bSHong Zhang   PetscFunctionBegin;
754ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
7626fbe8dcSKarl Rupp 
774ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
784ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
794ae0847bSHong Zhang   PetscFunctionReturn(0);
804ae0847bSHong Zhang }
814ae0847bSHong Zhang 
824ae0847bSHong Zhang #undef __FUNCT__
83a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
84a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
85598bc09dSHong Zhang {
86598bc09dSHong Zhang   PetscErrorCode ierr;
874ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
88598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
89598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
90598bc09dSHong Zhang   PetscInt       *adi=ad->i,*adj,*aoi=ao->i,*aoj;
91598bc09dSHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
92598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
93598bc09dSHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
94598bc09dSHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
95598bc09dSHong Zhang   PetscInt       cm   =C->rmap->n,anz,pnz;
96598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
9765a57a32SSean Farley   PetscInt       *api,*apj,*apJ,i,j,k,row;
9829825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
99598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
100*9467f45fSHong Zhang   MPI_Comm       comm;
101*9467f45fSHong Zhang   PetscMPIInt    size;
102598bc09dSHong Zhang 
103598bc09dSHong Zhang   PetscFunctionBegin;
104*9467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
105*9467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
106*9467f45fSHong Zhang 
107a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
108598bc09dSHong Zhang   /*-----------------------------------------------------*/
109a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
110b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
111a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
112598bc09dSHong Zhang 
113598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
114598bc09dSHong Zhang   /*----------------------------------------------------------*/
115598bc09dSHong Zhang   /* get data from symbolic products */
116a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
117598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
118*9467f45fSHong Zhang   if (size >1) {
119*9467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
120598bc09dSHong Zhang     pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
121*9467f45fSHong Zhang   }
122598bc09dSHong Zhang 
123598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
124598bc09dSHong Zhang   apa = ptap->apa;
125598bc09dSHong Zhang 
12629825010SHong Zhang   api = ptap->api;
12729825010SHong Zhang   apj = ptap->apj;
128598bc09dSHong Zhang   for (i=0; i<cm; i++) {
129598bc09dSHong Zhang     /* diagonal portion of A */
130598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
131598bc09dSHong Zhang     adj = ad->j + adi[i];
132598bc09dSHong Zhang     ada = ad->a + adi[i];
133598bc09dSHong Zhang     for (j=0; j<anz; j++) {
134598bc09dSHong Zhang       row = adj[j];
135598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
136598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
137598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
138598bc09dSHong Zhang 
139598bc09dSHong Zhang       /* perform dense axpy */
140598bc09dSHong Zhang       valtmp = ada[j];
141598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
142598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
143598bc09dSHong Zhang       }
144598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
145598bc09dSHong Zhang     }
146598bc09dSHong Zhang 
147598bc09dSHong Zhang     /* off-diagonal portion of A */
148598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
149598bc09dSHong Zhang     aoj = ao->j + aoi[i];
150598bc09dSHong Zhang     aoa = ao->a + aoi[i];
151598bc09dSHong Zhang     for (j=0; j<anz; j++) {
152598bc09dSHong Zhang       row = aoj[j];
153598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
154598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
155598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
156598bc09dSHong Zhang 
157598bc09dSHong Zhang       /* perform dense axpy */
158598bc09dSHong Zhang       valtmp = aoa[j];
159598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
160598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
161598bc09dSHong Zhang       }
162598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
163598bc09dSHong Zhang     }
164598bc09dSHong Zhang 
165598bc09dSHong Zhang     /* set values in C */
166598bc09dSHong Zhang     apJ  = apj + api[i];
167598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
168598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
169598bc09dSHong Zhang 
170598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
171598bc09dSHong Zhang     ca = coa + co->i[i];
172598bc09dSHong Zhang     k  = 0;
173598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
174598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
175598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
176598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
177598bc09dSHong Zhang       k++;
178598bc09dSHong Zhang     }
179598bc09dSHong Zhang 
180598bc09dSHong Zhang     /* diagonal part of C */
181598bc09dSHong Zhang     ca = cda + cd->i[i];
182598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
183598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
184598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
185598bc09dSHong Zhang       k++;
186598bc09dSHong Zhang     }
187598bc09dSHong Zhang 
188598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
189598bc09dSHong Zhang     ca = coa + co->i[i];
190598bc09dSHong Zhang     for (; k0<conz; k0++) {
191598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
192598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
193598bc09dSHong Zhang       k++;
194598bc09dSHong Zhang     }
195598bc09dSHong Zhang   }
196598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
197598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
198598bc09dSHong Zhang   PetscFunctionReturn(0);
199598bc09dSHong Zhang }
200598bc09dSHong Zhang 
201598bc09dSHong Zhang #undef __FUNCT__
2020fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
2030fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
204598bc09dSHong Zhang {
205598bc09dSHong Zhang   PetscErrorCode     ierr;
206ce94432eSBarry Smith   MPI_Comm           comm;
207*9467f45fSHong Zhang   PetscMPIInt        size;
208bfcd1a73SHong Zhang   Mat                Cmpi;
209598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
2100298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2114ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
212bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2134ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2144ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
215d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
216bfcd1a73SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
217598bc09dSHong Zhang   PetscBT            lnkbt;
218598bc09dSHong Zhang   PetscScalar        *apa;
219bfcd1a73SHong Zhang   PetscReal          afill;
220f84c536bSHong Zhang   PetscInt           nlnk_max,armax,prmax;
221598bc09dSHong Zhang 
222598bc09dSHong Zhang   PetscFunctionBegin;
223ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
224*9467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
225*9467f45fSHong Zhang 
226a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
227cf1a0672SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
228a1a4e84aSHong Zhang   }
229152983bfSHong Zhang 
230598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
231b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
232598bc09dSHong Zhang 
233598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
234b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
23519f0e57aSHong Zhang 
236598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
237a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
238598bc09dSHong Zhang 
239a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
240598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
241*9467f45fSHong Zhang   if (size > 1) {
242*9467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
243598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
244*9467f45fSHong Zhang   }
245598bc09dSHong Zhang 
246598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
247598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
248785e854fSJed Brown   ierr      = PetscMalloc1((am+2),&api);CHKERRQ(ierr);
249a1a4e84aSHong Zhang   ptap->api = api;
250598bc09dSHong Zhang   api[0]    = 0;
251598bc09dSHong Zhang 
252598bc09dSHong Zhang   /* create and initialize a linked list */
253f84c536bSHong Zhang   armax    = ad->rmax+ao->rmax;
254*9467f45fSHong Zhang   if (size >1) {
255f84c536bSHong Zhang     prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
256*9467f45fSHong Zhang   } else {
257*9467f45fSHong Zhang     prmax = p_loc->rmax;
258*9467f45fSHong Zhang   }
259f84c536bSHong Zhang   nlnk_max = armax*prmax;
260f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2610d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
262598bc09dSHong Zhang 
263bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
264bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
2652205254eSKarl Rupp 
266598bc09dSHong Zhang   current_space = free_space;
267598bc09dSHong Zhang 
268598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
269598bc09dSHong Zhang   for (i=0; i<am; i++) {
270598bc09dSHong Zhang     /* diagonal portion of A */
271598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
272598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
273598bc09dSHong Zhang       row  = *adj++;
274598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
275598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
276598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2771fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
278598bc09dSHong Zhang     }
279598bc09dSHong Zhang     /* off-diagonal portion of A */
280598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
281598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
282598bc09dSHong Zhang       row  = *aoj++;
283598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
284598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2851fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
286598bc09dSHong Zhang     }
287598bc09dSHong Zhang 
288d14fa243SHong Zhang     apnz     = lnk[0];
289598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
290598bc09dSHong Zhang 
291598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
292598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
293598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
294598bc09dSHong Zhang       nspacedouble++;
295598bc09dSHong Zhang     }
296598bc09dSHong Zhang 
297598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2981fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
299598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
3002205254eSKarl Rupp 
301598bc09dSHong Zhang     current_space->array           += apnz;
302598bc09dSHong Zhang     current_space->local_used      += apnz;
303598bc09dSHong Zhang     current_space->local_remaining -= apnz;
304598bc09dSHong Zhang   }
305598bc09dSHong Zhang 
306598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
307598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
308785e854fSJed Brown   ierr = PetscMalloc1((api[am]+1),&ptap->apj);CHKERRQ(ierr);
309a1a4e84aSHong Zhang   apj  = ptap->apj;
310a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
311598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
312598bc09dSHong Zhang 
313f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
3141795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
3152205254eSKarl Rupp 
316f84c536bSHong Zhang   ptap->apa = apa;
317f84c536bSHong Zhang 
318bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
319598bc09dSHong Zhang   /*----------------------------------------------------*/
320bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
321bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
32233d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
323a2f3521dSMark F. Adams 
324bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
325bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
326598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
327598bc09dSHong Zhang   for (i=0; i<am; i++) {
328598bc09dSHong Zhang     row  = i + rstart;
329598bc09dSHong Zhang     apnz = api[i+1] - api[i];
330bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
331598bc09dSHong Zhang     apj += apnz;
332598bc09dSHong Zhang   }
333bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
334bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
335598bc09dSHong Zhang 
336bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
337bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
338bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
339bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
340598bc09dSHong Zhang 
341bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
342bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
343598bc09dSHong Zhang   c->ptap = ptap;
344598bc09dSHong Zhang 
345bfcd1a73SHong Zhang   *C = Cmpi;
346bfcd1a73SHong Zhang 
347bfcd1a73SHong Zhang   /* set MatInfo */
348118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
349bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
350bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
351bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
352bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
353bfcd1a73SHong Zhang 
354bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
355bfcd1a73SHong Zhang   if (api[am]) {
35657622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
35757622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
358bfcd1a73SHong Zhang   } else {
359bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
360bfcd1a73SHong Zhang   }
361bfcd1a73SHong Zhang #endif
362598bc09dSHong Zhang   PetscFunctionReturn(0);
363598bc09dSHong Zhang }
364598bc09dSHong Zhang 
3659123193aSHong Zhang #undef __FUNCT__
3669123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3679123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3689123193aSHong Zhang {
3699123193aSHong Zhang   PetscErrorCode ierr;
3709123193aSHong Zhang 
3719123193aSHong Zhang   PetscFunctionBegin;
3729123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3733ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3749123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3753ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3769123193aSHong Zhang   }
3773ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3789123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3793ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3809123193aSHong Zhang   PetscFunctionReturn(0);
3819123193aSHong Zhang }
3829123193aSHong Zhang 
3834324174eSBarry Smith typedef struct {
3844324174eSBarry Smith   Mat         workB;
3854324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3864324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3874324174eSBarry Smith } MPIAIJ_MPIDense;
3884324174eSBarry Smith 
3897af9e4fdSHong Zhang #undef __FUNCT__
39096b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
39196b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3924324174eSBarry Smith {
3934324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3944324174eSBarry Smith   PetscErrorCode  ierr;
3954324174eSBarry Smith 
3964324174eSBarry Smith   PetscFunctionBegin;
3976bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3984324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3994324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
4004324174eSBarry Smith   PetscFunctionReturn(0);
4014324174eSBarry Smith }
4024324174eSBarry Smith 
4039123193aSHong Zhang #undef __FUNCT__
4049123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4059123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4069123193aSHong Zhang {
4079123193aSHong Zhang   PetscErrorCode         ierr;
408f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
409d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
410bf0cc555SLisandro Dalcin   PetscContainer         container;
4114324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4124324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4134324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4144324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
415d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4169123193aSHong Zhang 
4179123193aSHong Zhang   PetscFunctionBegin;
418ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
419d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
42033d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
421cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4220298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
423cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
424cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4252205254eSKarl Rupp 
426d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
427f32d5d43SBarry Smith 
428b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
429f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4300298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4314324174eSBarry Smith   /* Create work arrays needed */
432dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
433dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
434dcca6d9dSJed Brown                       from->n,&contents->rwaits,
435dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4364324174eSBarry Smith 
437ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
438bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
43996b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
440bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
441bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4429123193aSHong Zhang   PetscFunctionReturn(0);
4439123193aSHong Zhang }
4449123193aSHong Zhang 
4457af9e4fdSHong Zhang #undef __FUNCT__
4467af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
447f32d5d43SBarry Smith /*
4482636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4492636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
450f32d5d43SBarry Smith */
4514324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
452f32d5d43SBarry Smith {
453f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
454f32d5d43SBarry Smith   PetscErrorCode         ierr;
455f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
456f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
457f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
458f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
459f32d5d43SBarry Smith   PetscInt               i,j,k;
460aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
461aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
462f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
463ce94432eSBarry Smith   MPI_Comm               comm;
464d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
465f32d5d43SBarry Smith   MPI_Status             status;
4664324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
467bf0cc555SLisandro Dalcin   PetscContainer         container;
4684324174eSBarry Smith   Mat                    workB;
469f32d5d43SBarry Smith 
470f32d5d43SBarry Smith   PetscFunctionBegin;
471ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
472bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
47329825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
474bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4754324174eSBarry Smith 
4764324174eSBarry Smith   workB = *outworkB = contents->workB;
477cf1a0672SHong 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);
478f32d5d43SBarry Smith   sindices = to->indices;
479f32d5d43SBarry Smith   sstarts  = to->starts;
480f32d5d43SBarry Smith   sprocs   = to->procs;
4814324174eSBarry Smith   swaits   = contents->swaits;
4824324174eSBarry Smith   svalues  = contents->svalues;
483f32d5d43SBarry Smith 
484f32d5d43SBarry Smith   rindices = from->indices;
485f32d5d43SBarry Smith   rstarts  = from->starts;
486f32d5d43SBarry Smith   rprocs   = from->procs;
4874324174eSBarry Smith   rwaits   = contents->rwaits;
4884324174eSBarry Smith   rvalues  = contents->rvalues;
489f32d5d43SBarry Smith 
4908c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
4918c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
492f32d5d43SBarry Smith 
493f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
494f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
495f32d5d43SBarry Smith   }
496f32d5d43SBarry Smith 
497f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
498f32d5d43SBarry Smith     /* pack a message at a time */
499f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
500f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
501f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5022636ff26SBarry Smith       }
5032636ff26SBarry Smith     }
504f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
505f32d5d43SBarry Smith   }
5062636ff26SBarry Smith 
507f32d5d43SBarry Smith   nrecvs = from->n;
508f32d5d43SBarry Smith   while (nrecvs) {
509f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
510f32d5d43SBarry Smith     nrecvs--;
511f32d5d43SBarry Smith     /* unpack a message at a time */
512f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
513f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
514f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5152636ff26SBarry Smith       }
5162636ff26SBarry Smith     }
517f32d5d43SBarry Smith   }
518cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
519f32d5d43SBarry Smith 
5208c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5218c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
522f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
523f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
524f32d5d43SBarry Smith   PetscFunctionReturn(0);
525f32d5d43SBarry Smith }
526f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
527f32d5d43SBarry Smith 
5289123193aSHong Zhang #undef __FUNCT__
5299123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5309123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5319123193aSHong Zhang {
5329123193aSHong Zhang   PetscErrorCode ierr;
533f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
534f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
535f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
536f32d5d43SBarry Smith   Mat            workB;
5379123193aSHong Zhang 
5389123193aSHong Zhang   PetscFunctionBegin;
539f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
540f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
541f32d5d43SBarry Smith 
542f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5434324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
544f32d5d43SBarry Smith 
545f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
546f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5479123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5489123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5499123193aSHong Zhang   PetscFunctionReturn(0);
5509123193aSHong Zhang }
551cf1a0672SHong Zhang 
5521634b032SHong Zhang #undef __FUNCT__
553cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
554cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
5551634b032SHong Zhang {
5561634b032SHong Zhang   PetscErrorCode ierr;
557cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
558cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
559cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
560cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
561cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
562cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
563cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
56429825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
56529825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
566cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
56729825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
568cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
56929825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
57029825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
571a7c7454dSHong Zhang   MPI_Comm       comm;
572a7c7454dSHong Zhang   PetscMPIInt    size;
5731634b032SHong Zhang 
5741634b032SHong Zhang   PetscFunctionBegin;
575a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
576a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
577a7c7454dSHong Zhang 
578cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
579cf1a0672SHong Zhang   /*-----------------------------------------------------*/
580cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
581b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
582cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
583cf1a0672SHong Zhang 
584cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
585cf1a0672SHong Zhang   /*----------------------------------------------------------*/
586cf1a0672SHong Zhang   /* get data from symbolic products */
587cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
588cf1a0672SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
589a7c7454dSHong Zhang   if (size >1) {
590a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
591cf1a0672SHong Zhang     pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
592a7c7454dSHong Zhang   }
593cf1a0672SHong Zhang 
59429825010SHong Zhang   api = ptap->api;
59529825010SHong Zhang   apj = ptap->apj;
596cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
59729825010SHong Zhang     apJ = apj + api[i];
59829825010SHong Zhang 
599cf1a0672SHong Zhang     /* diagonal portion of A */
600cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
601cf1a0672SHong Zhang     adj = ad->j + adi[i];
602cf1a0672SHong Zhang     ada = ad->a + adi[i];
603cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
604cf1a0672SHong Zhang       row = adj[j];
605cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
606cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
607cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
60829825010SHong Zhang       /* perform sparse axpy */
609cf1a0672SHong Zhang       valtmp = ada[j];
61029825010SHong Zhang       nextp  = 0;
61129825010SHong Zhang       for (k=0; nextp<pnz; k++) {
61229825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
61329825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
61429825010SHong Zhang         }
6151634b032SHong Zhang       }
616cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
617cf1a0672SHong Zhang     }
6181634b032SHong Zhang 
619cf1a0672SHong Zhang     /* off-diagonal portion of A */
620cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
621cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
622cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
623cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
624cf1a0672SHong Zhang       row = aoj[j];
625cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
626cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
627cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
62829825010SHong Zhang       /* perform sparse axpy */
629cf1a0672SHong Zhang       valtmp = aoa[j];
63029825010SHong Zhang       nextp  = 0;
63129825010SHong Zhang       for (k=0; nextp<pnz; k++) {
63229825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
63329825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
63429825010SHong Zhang         }
635cf1a0672SHong Zhang       }
636cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
637cf1a0672SHong Zhang     }
6381634b032SHong Zhang 
639cf1a0672SHong Zhang     /* set values in C */
640cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
641cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6421634b032SHong Zhang 
643cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
644cf1a0672SHong Zhang     ca = coa + co->i[i];
645cf1a0672SHong Zhang     k  = 0;
646cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
647cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
64829825010SHong Zhang       ca[k0]        = apa_sparse[k];
64929825010SHong Zhang       apa_sparse[k] = 0.0;
650cf1a0672SHong Zhang       k++;
651cf1a0672SHong Zhang     }
6521634b032SHong Zhang 
653cf1a0672SHong Zhang     /* diagonal part of C */
654cf1a0672SHong Zhang     ca = cda + cd->i[i];
655cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
65629825010SHong Zhang       ca[k1]        = apa_sparse[k];
65729825010SHong Zhang       apa_sparse[k] = 0.0;
658cf1a0672SHong Zhang       k++;
659cf1a0672SHong Zhang     }
660cf1a0672SHong Zhang 
661cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
662cf1a0672SHong Zhang     ca = coa + co->i[i];
663cf1a0672SHong Zhang     for (; k0<conz; k0++) {
66429825010SHong Zhang       ca[k0]        = apa_sparse[k];
66529825010SHong Zhang       apa_sparse[k] = 0.0;
666cf1a0672SHong Zhang       k++;
667cf1a0672SHong Zhang     }
668cf1a0672SHong Zhang   }
669cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
670cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
671cf1a0672SHong Zhang   PetscFunctionReturn(0);
672cf1a0672SHong Zhang }
673cf1a0672SHong Zhang 
6740fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
675cf1a0672SHong Zhang #undef __FUNCT__
676b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
677b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
678cf1a0672SHong Zhang {
679cf1a0672SHong Zhang   PetscErrorCode     ierr;
680ce94432eSBarry Smith   MPI_Comm           comm;
681a7c7454dSHong Zhang   PetscMPIInt        size;
682cf1a0672SHong Zhang   Mat                Cmpi;
683cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6840298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
685cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
686cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
687cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
688cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
689f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
690cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
691cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
692cf1a0672SHong Zhang   PetscReal          afill;
693cf1a0672SHong Zhang   PetscScalar        *apa;
694cf1a0672SHong Zhang 
695cf1a0672SHong Zhang   PetscFunctionBegin;
696ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
697a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
698a7c7454dSHong Zhang 
699cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
700b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
701cf1a0672SHong Zhang 
702cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
703b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
70419f0e57aSHong Zhang 
705cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
706cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
707cf1a0672SHong Zhang 
708cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
709cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
710a7c7454dSHong Zhang   if (size > 1) {
711a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
712a7c7454dSHong Zhang     pi_oth = p_oth->i;
713a7c7454dSHong Zhang     pj_oth = p_oth->j;
714a7c7454dSHong Zhang   }
715cf1a0672SHong Zhang 
716cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
717cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
718785e854fSJed Brown   ierr      = PetscMalloc1((am+2),&api);CHKERRQ(ierr);
719cf1a0672SHong Zhang   ptap->api = api;
720cf1a0672SHong Zhang   api[0]    = 0;
721cf1a0672SHong Zhang 
722cf1a0672SHong Zhang   /* create and initialize a linked list */
723cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
724a7c7454dSHong Zhang   if (size >1) {
725cf1a0672SHong Zhang     prmax = PetscMax(p_loc->rmax,p_oth->rmax);
726a7c7454dSHong Zhang   } else {
727a7c7454dSHong Zhang     prmax = p_loc->rmax;
728a7c7454dSHong Zhang   }
729cf1a0672SHong Zhang   nlnk_max = armax*prmax;
730cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
731f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
732cf1a0672SHong Zhang 
733cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
734cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
7352205254eSKarl Rupp 
736cf1a0672SHong Zhang   current_space = free_space;
737cf1a0672SHong Zhang 
738cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
739cf1a0672SHong Zhang   for (i=0; i<am; i++) {
740cf1a0672SHong Zhang     /* diagonal portion of A */
741cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
742cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
743cf1a0672SHong Zhang       row  = *adj++;
744cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
745cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
746cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
747f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
748cf1a0672SHong Zhang     }
749cf1a0672SHong Zhang     /* off-diagonal portion of A */
750cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
751cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
752cf1a0672SHong Zhang       row  = *aoj++;
753cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
754cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
755f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
756cf1a0672SHong Zhang     }
757cf1a0672SHong Zhang 
758f84c536bSHong Zhang     apnz     = *lnk;
759cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
760cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
761cf1a0672SHong Zhang 
762cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
763cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
764cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
765cf1a0672SHong Zhang       nspacedouble++;
766cf1a0672SHong Zhang     }
767cf1a0672SHong Zhang 
768cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
769f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
770cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7712205254eSKarl Rupp 
772cf1a0672SHong Zhang     current_space->array           += apnz;
773cf1a0672SHong Zhang     current_space->local_used      += apnz;
774cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
775cf1a0672SHong Zhang   }
776cf1a0672SHong Zhang 
777cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
778cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
779785e854fSJed Brown   ierr = PetscMalloc1((api[am]+1),&ptap->apj);CHKERRQ(ierr);
780cf1a0672SHong Zhang   apj  = ptap->apj;
781cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
782f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
783cf1a0672SHong Zhang 
784cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
785cf1a0672SHong Zhang   /*----------------------------------------------------*/
786cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
787cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
78833d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
789cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
790cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
791cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
792cf1a0672SHong Zhang 
79329825010SHong Zhang   /* malloc apa for assembly Cmpi */
7941795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
7952205254eSKarl Rupp 
79629825010SHong Zhang   ptap->apa = apa;
797cf1a0672SHong Zhang   for (i=0; i<am; i++) {
798cf1a0672SHong Zhang     row  = i + rstart;
799cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
800cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
801cf1a0672SHong Zhang     apj += apnz;
802cf1a0672SHong Zhang   }
803cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
804cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
805cf1a0672SHong Zhang 
806cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
807cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
808cf1a0672SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
809cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
810cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
811cf1a0672SHong Zhang 
812cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
813cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
814cf1a0672SHong Zhang   c->ptap = ptap;
815cf1a0672SHong Zhang 
816cf1a0672SHong Zhang   *C = Cmpi;
817cf1a0672SHong Zhang 
818cf1a0672SHong Zhang   /* set MatInfo */
819118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
820cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
821cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
822cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
823cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
824cf1a0672SHong Zhang 
825cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
826cf1a0672SHong Zhang   if (api[am]) {
82757622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
82857622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
829cf1a0672SHong Zhang   } else {
830cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
831cf1a0672SHong Zhang   }
832cf1a0672SHong Zhang #endif
8331634b032SHong Zhang   PetscFunctionReturn(0);
8341634b032SHong Zhang }
835f96d37fcSHong Zhang 
836f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
837f96d37fcSHong Zhang #undef __FUNCT__
838f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
839c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
840f96d37fcSHong Zhang {
841f96d37fcSHong Zhang   PetscErrorCode ierr;
842c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
843c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
844f96d37fcSHong Zhang 
845f96d37fcSHong Zhang   PetscFunctionBegin;
846c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
847c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
848c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
849c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
850c216dbf3SHong Zhang 
851c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
852c216dbf3SHong Zhang     switch (alg) {
853c216dbf3SHong Zhang     case 1:
8542bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
855c216dbf3SHong Zhang       break;
856c216dbf3SHong Zhang     case 2:
857275476c6SMatthew G. Knepley     {
858ab1b013aSHong Zhang       Mat         Pt;
859ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
860ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
861ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
862ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
863ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
864ab1b013aSHong Zhang       ptap     = c->ptap;
865ab1b013aSHong Zhang       ptap->Pt = Pt;
866c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
867c216dbf3SHong Zhang       PetscFunctionReturn(0);
868275476c6SMatthew G. Knepley     }
869c216dbf3SHong Zhang       break;
870c216dbf3SHong Zhang     default:
8716da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
872c216dbf3SHong Zhang       break;
873c216dbf3SHong Zhang     }
874c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
875c216dbf3SHong Zhang   }
876c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
877c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
878c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
879ab1b013aSHong Zhang   PetscFunctionReturn(0);
880ab1b013aSHong Zhang }
881ab1b013aSHong Zhang 
882c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
883c216dbf3SHong Zhang #undef __FUNCT__
884c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
885c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
886c216dbf3SHong Zhang {
887c216dbf3SHong Zhang   PetscErrorCode ierr;
8882bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
8892bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
8902bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
891c216dbf3SHong Zhang 
892c216dbf3SHong Zhang   PetscFunctionBegin;
893c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
894c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
895f96d37fcSHong Zhang   PetscFunctionReturn(0);
896f96d37fcSHong Zhang }
897f96d37fcSHong Zhang 
8986da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
899f96d37fcSHong Zhang #undef __FUNCT__
9006da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9016da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
902f96d37fcSHong Zhang {
903c5af039cSHong Zhang   PetscErrorCode      ierr;
904c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
905497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
906c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
907c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
908d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
909497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
910e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
911c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
912ce94432eSBarry Smith   MPI_Comm            comm;
913c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
914c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
915c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
916c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
917c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
918c5af039cSHong Zhang   MPI_Status          *status;
919c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
920d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
921a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
922c5af039cSHong Zhang   Mat                 A_loc;
923c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
924f96d37fcSHong Zhang 
925f96d37fcSHong Zhang   PetscFunctionBegin;
926ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
927c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
928c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
929c5af039cSHong Zhang 
930c5af039cSHong Zhang   ptap  = c->ptap;
931c5af039cSHong Zhang   merge = ptap->merge;
932c5af039cSHong Zhang 
933c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
934c5af039cSHong Zhang   /*--------------------------------------------------------------*/
935c5af039cSHong Zhang   /* get data from symbolic products */
936c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
9371795a4d1SJed Brown   ierr = PetscCalloc1((coi[pon]+1),&coa);CHKERRQ(ierr);
938c5af039cSHong Zhang 
939c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
940c5af039cSHong Zhang   owners = merge->rowmap->range;
9411795a4d1SJed Brown   ierr   = PetscCalloc1((bi[cm]+1),&ba);CHKERRQ(ierr);
942c5af039cSHong Zhang 
943c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
944c5af039cSHong Zhang   A_loc = ptap->A_loc;
945c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
946c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
947d6ab1dc0SHong Zhang   ai    = a_loc->i;
948d6ab1dc0SHong Zhang   aj    = a_loc->j;
949c5af039cSHong Zhang 
9501795a4d1SJed Brown   ierr = PetscCalloc1((A->cmap->N),&aval);CHKERRQ(ierr); /* non-scalable!!! */
951c5af039cSHong Zhang 
952c5af039cSHong Zhang   for (i=0; i<am; i++) {
953e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
954d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
955d6ab1dc0SHong Zhang     adj = aj + ai[i];
956d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
957c5af039cSHong Zhang     for (j=0; j<anz; j++) {
958e745005bSHong Zhang       aval[adj[j]] = ada[j];
959c5af039cSHong Zhang     }
960c5af039cSHong Zhang 
961c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
962c5af039cSHong Zhang     /*--------------------------------------------------------------*/
963c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
964c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
965c5af039cSHong Zhang     poJ = po->j + po->i[i];
966c5af039cSHong Zhang     pA  = po->a + po->i[i];
967c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
968c5af039cSHong Zhang       row = poJ[j];
969c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
970c5af039cSHong Zhang       cj  = coj + coi[row];
971c5af039cSHong Zhang       ca  = coa + coi[row];
972c5af039cSHong Zhang       /* perform dense axpy */
973c5af039cSHong Zhang       valtmp = pA[j];
974c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
975e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
976c5af039cSHong Zhang       }
977c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
978c5af039cSHong Zhang     }
979c5af039cSHong Zhang 
980c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
981c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
982c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
983c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
984c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
985c5af039cSHong Zhang       row = pdJ[j];
986c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
987c5af039cSHong Zhang       cj  = bj + bi[row];
988c5af039cSHong Zhang       ca  = ba + bi[row];
989c5af039cSHong Zhang       /* perform dense axpy */
990c5af039cSHong Zhang       valtmp = pA[j];
991c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
992e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
993c5af039cSHong Zhang       }
994c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
995c5af039cSHong Zhang     }
996c5af039cSHong Zhang 
997d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
998d6ab1dc0SHong Zhang     aJ = aj + ai[i];
999e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1000c5af039cSHong Zhang   }
1001c5af039cSHong Zhang 
1002c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1003c5af039cSHong Zhang   /*------------------------------------*/
1004c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1005c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1006c5af039cSHong Zhang   len_s  = merge->len_s;
1007c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1008c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1009c5af039cSHong Zhang 
1010dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1011c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1012c5af039cSHong Zhang     if (!len_s[proc]) continue;
1013c5af039cSHong Zhang     i    = merge->owners_co[proc];
1014c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1015c5af039cSHong Zhang     k++;
1016c5af039cSHong Zhang   }
1017c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1018c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1019c5af039cSHong Zhang 
1020c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1021c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1022c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1023c5af039cSHong Zhang 
1024c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1025c5af039cSHong Zhang   /*----------------------------------------------------*/
1026dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1027c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1028c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1029c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1030c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1031c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1032c5af039cSHong Zhang   }
1033c5af039cSHong Zhang 
1034c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1035c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1036c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1037c5af039cSHong Zhang     ba_i = ba + bi[i];
1038c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1039c5af039cSHong Zhang     /* add received vals into ba */
1040c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1041c5af039cSHong Zhang       /* i-th row */
1042c5af039cSHong Zhang       if (i == *nextrow[k]) {
1043c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1044c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1045c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1046c5af039cSHong Zhang         nextcj = 0;
1047c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1048c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1049c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1050c5af039cSHong Zhang           }
1051c5af039cSHong Zhang         }
1052c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1053c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1054c5af039cSHong Zhang       }
1055c5af039cSHong Zhang     }
1056c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1057c5af039cSHong Zhang   }
1058c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1059c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1060c5af039cSHong Zhang 
1061c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1062c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1063c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1064c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1065e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1066f96d37fcSHong Zhang   PetscFunctionReturn(0);
1067f96d37fcSHong Zhang }
1068f96d37fcSHong Zhang 
1069c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1070f96d37fcSHong Zhang #undef __FUNCT__
10712bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
10722bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1073f96d37fcSHong Zhang {
1074f96d37fcSHong Zhang   PetscErrorCode      ierr;
10754e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1076f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10770298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1078f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1079f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1080f96d37fcSHong Zhang   PetscInt            nnz;
10814e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1082497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1083f96d37fcSHong Zhang   PetscBT             lnkbt;
1084ce94432eSBarry Smith   MPI_Comm            comm;
1085f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1086f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1087f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1088f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1089f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1090f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1091f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1092f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1093d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1094f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1095c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1096f96d37fcSHong Zhang   PetscScalar         *vals;
10974e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1098f96d37fcSHong Zhang 
1099f96d37fcSHong Zhang   PetscFunctionBegin;
1100ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1101c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1102c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1103c5af039cSHong 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);
1104c5af039cSHong Zhang   }
1105c5af039cSHong Zhang 
1106f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1107f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1108f96d37fcSHong Zhang 
1109f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1110b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1111f96d37fcSHong Zhang 
1112f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1113f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11142205254eSKarl Rupp 
1115c5af039cSHong Zhang   ptap->A_loc = A_loc;
11162205254eSKarl Rupp 
11171c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1118d6ab1dc0SHong Zhang   ai    = a_loc->i;
1119d6ab1dc0SHong Zhang   aj    = a_loc->j;
1120f96d37fcSHong Zhang 
1121f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1122f96d37fcSHong Zhang   /*----------------------------------------------------*/
11234e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11244e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11254e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11264e938277SHong Zhang 
11274e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11284e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11294e938277SHong Zhang   poti = pot->i; potj = pot->j;
1130f96d37fcSHong Zhang 
1131f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11322205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1133785e854fSJed Brown   ierr   = PetscMalloc1((pon+1),&coi);CHKERRQ(ierr);
1134f96d37fcSHong Zhang   coi[0] = 0;
1135f96d37fcSHong Zhang 
1136f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1137d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1138a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1139f96d37fcSHong Zhang   current_space = free_space;
114019f0e57aSHong Zhang 
1141f96d37fcSHong Zhang   /* create and initialize a linked list */
11424e938277SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1143c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
11444e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
11454e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1146f96d37fcSHong Zhang 
1147f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1148f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1149f96d37fcSHong Zhang     ptJ = potj + poti[i];
1150f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1151f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1152d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1153d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1154f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1155d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1156f96d37fcSHong Zhang     }
11574e938277SHong Zhang     nnz = lnk[0];
1158f96d37fcSHong Zhang 
1159f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1160f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1161f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1162f96d37fcSHong Zhang       nspacedouble++;
1163f96d37fcSHong Zhang     }
1164f96d37fcSHong Zhang 
1165f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11664e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11672205254eSKarl Rupp 
1168f96d37fcSHong Zhang     current_space->array           += nnz;
1169f96d37fcSHong Zhang     current_space->local_used      += nnz;
1170f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11712205254eSKarl Rupp 
1172f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1173f96d37fcSHong Zhang   }
1174f96d37fcSHong Zhang 
1175785e854fSJed Brown   ierr = PetscMalloc1((coi[pon]+1),&coj);CHKERRQ(ierr);
1176f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11772205254eSKarl Rupp 
1178118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1179f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1180f96d37fcSHong Zhang 
1181f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1182f96d37fcSHong Zhang   /*----------------------------------------------*/
1183f96d37fcSHong Zhang   /* determine row ownership */
1184b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1185f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11862205254eSKarl Rupp 
1187f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1188f96d37fcSHong Zhang   merge->rowmap->bs = 1;
11892205254eSKarl Rupp 
1190f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1191f96d37fcSHong Zhang   owners = merge->rowmap->range;
1192f96d37fcSHong Zhang 
1193f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
11941795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1195785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
11962205254eSKarl Rupp 
1197f96d37fcSHong Zhang   len_s        = merge->len_s;
1198f96d37fcSHong Zhang   merge->nsend = 0;
1199f96d37fcSHong Zhang 
1200785e854fSJed Brown   ierr = PetscMalloc1((size+2),&owners_co);CHKERRQ(ierr);
1201f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1202f96d37fcSHong Zhang 
1203f96d37fcSHong Zhang   proc = 0;
1204f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1205f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1206f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1207f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1208f96d37fcSHong Zhang   }
1209f96d37fcSHong Zhang 
1210f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1211f96d37fcSHong Zhang   owners_co[0] = 0;
1212f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1213f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1214f96d37fcSHong Zhang     if (len_si[proc]) {
1215f96d37fcSHong Zhang       merge->nsend++;
1216f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1217f96d37fcSHong Zhang       len         += len_si[proc];
1218f96d37fcSHong Zhang     }
1219f96d37fcSHong Zhang   }
1220f96d37fcSHong Zhang 
1221f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12220298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1223f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1224f96d37fcSHong Zhang 
1225f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1226f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1227f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1228785e854fSJed Brown   ierr = PetscMalloc1((merge->nsend+1),&swaits);CHKERRQ(ierr);
1229f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1230f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1231f96d37fcSHong Zhang     i    = owners_co[proc];
1232f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1233f96d37fcSHong Zhang     k++;
1234f96d37fcSHong Zhang   }
1235f96d37fcSHong Zhang 
1236f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1237785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1238f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1239c280ed6aSJed Brown     PetscMPIInt icompleted;
1240c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1241f96d37fcSHong Zhang   }
1242f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1243f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1244f96d37fcSHong Zhang 
1245f96d37fcSHong Zhang   /* send and recv coi */
1246f96d37fcSHong Zhang   /*-------------------*/
1247f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1248f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1249785e854fSJed Brown   ierr   = PetscMalloc1((len+1),&buf_s);CHKERRQ(ierr);
1250f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1251f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1252f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1253f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1254f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1255f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1256f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1257f96d37fcSHong Zhang     */
1258f96d37fcSHong Zhang     /*-------------------------------------------*/
1259f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1260f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1261f96d37fcSHong Zhang     buf_si[0]   = nrows;
1262f96d37fcSHong Zhang     buf_si_i[0] = 0;
1263f96d37fcSHong Zhang     nrows       = 0;
1264f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1265f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1266f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1267f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1268f96d37fcSHong Zhang       nrows++;
1269f96d37fcSHong Zhang     }
1270f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1271f96d37fcSHong Zhang     k++;
1272f96d37fcSHong Zhang     buf_si += len_si[proc];
1273f96d37fcSHong Zhang   }
1274f96d37fcSHong Zhang   i = merge->nrecv;
1275f96d37fcSHong Zhang   while (i--) {
1276c280ed6aSJed Brown     PetscMPIInt icompleted;
1277c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1278f96d37fcSHong Zhang   }
1279f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1280f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1281f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1282f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1283f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1284f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1285f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1286f96d37fcSHong Zhang 
1287f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1288f96d37fcSHong Zhang   /*------------------------------------------*/
1289f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1290785e854fSJed Brown   ierr  = PetscMalloc1((pn+1),&bi);CHKERRQ(ierr);
1291f96d37fcSHong Zhang   bi[0] = 0;
1292f96d37fcSHong Zhang 
1293c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1294d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1295a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1296f96d37fcSHong Zhang   current_space = free_space;
1297f96d37fcSHong Zhang 
1298dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1299f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1300f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1301f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1302f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1303f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1304f96d37fcSHong Zhang   }
1305f96d37fcSHong Zhang 
13061c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1307f96d37fcSHong Zhang   rmax = 0;
1308f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1309f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1310f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1311f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1312f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1313f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1314d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1315d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1316f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1317d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1318f96d37fcSHong Zhang     }
13194e938277SHong Zhang 
1320f96d37fcSHong Zhang     /* add received col data into lnk */
1321f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1322f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1323f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1324f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13254e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1326f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1327f96d37fcSHong Zhang       }
1328f96d37fcSHong Zhang     }
13294e938277SHong Zhang     nnz = lnk[0];
1330f96d37fcSHong Zhang 
1331f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1332f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1333f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1334f96d37fcSHong Zhang       nspacedouble++;
1335f96d37fcSHong Zhang     }
1336f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13374e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1338f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13392205254eSKarl Rupp 
1340f96d37fcSHong Zhang     current_space->array           += nnz;
1341f96d37fcSHong Zhang     current_space->local_used      += nnz;
1342f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13432205254eSKarl Rupp 
1344f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1345f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1346f96d37fcSHong Zhang   }
1347f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1348f96d37fcSHong Zhang 
1349785e854fSJed Brown   ierr = PetscMalloc1((bi[pn]+1),&bj);CHKERRQ(ierr);
1350f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13512205254eSKarl Rupp 
1352118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1353f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1354d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13554e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13564e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1357f96d37fcSHong Zhang 
13581c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13591c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13601795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1361f96d37fcSHong Zhang 
1362f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1363f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
136433d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1365f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1366f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1367f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1368f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1369f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1370f96d37fcSHong Zhang     row  = i + rstart;
1371f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1372f96d37fcSHong Zhang     Jptr = bj + bi[i];
1373f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1374f96d37fcSHong Zhang   }
1375f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1376f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1377f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1378f96d37fcSHong Zhang 
1379f96d37fcSHong Zhang   merge->bi        = bi;
1380f96d37fcSHong Zhang   merge->bj        = bj;
1381f96d37fcSHong Zhang   merge->coi       = coi;
1382f96d37fcSHong Zhang   merge->coj       = coj;
1383f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1384f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1385f96d37fcSHong Zhang   merge->owners_co = owners_co;
1386f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1387f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1388f96d37fcSHong Zhang 
13896da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1390c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1391f96d37fcSHong Zhang 
1392f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1393f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1394f96d37fcSHong Zhang   c->ptap     = ptap;
13950298fd71SBarry Smith   ptap->api   = NULL;
13960298fd71SBarry Smith   ptap->apj   = NULL;
1397f96d37fcSHong Zhang   ptap->merge = merge;
1398d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1399d6ab1dc0SHong Zhang 
1400d6ab1dc0SHong Zhang   *C = Cmpi;
1401d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1402d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
140357622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
140457622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1405d6ab1dc0SHong Zhang   } else {
1406d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1407d6ab1dc0SHong Zhang   }
1408d6ab1dc0SHong Zhang #endif
1409d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1410d6ab1dc0SHong Zhang }
1411d6ab1dc0SHong Zhang 
1412d6ab1dc0SHong Zhang #undef __FUNCT__
14136da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14146da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1415d6ab1dc0SHong Zhang {
1416d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1417d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1418d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1419d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1420d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1421e745005bSHong Zhang   PetscInt            *adj;
1422e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1423e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1424d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1425ce94432eSBarry Smith   MPI_Comm            comm;
1426d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1427d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1428d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1429d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1430d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1431d6ab1dc0SHong Zhang   MPI_Status          *status;
1432d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1433d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1434a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1435d6ab1dc0SHong Zhang   Mat                 A_loc;
1436d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1437d6ab1dc0SHong Zhang 
1438d6ab1dc0SHong Zhang   PetscFunctionBegin;
1439ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1440d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1441d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1442d6ab1dc0SHong Zhang 
1443d6ab1dc0SHong Zhang   ptap  = c->ptap;
1444d6ab1dc0SHong Zhang   merge = ptap->merge;
1445d6ab1dc0SHong Zhang 
1446e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1447e745005bSHong Zhang   /*------------------------------------------*/
1448d6ab1dc0SHong Zhang   /* get data from symbolic products */
1449d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
14501795a4d1SJed Brown   ierr   = PetscCalloc1((coi[pon]+1),&coa);CHKERRQ(ierr);
1451d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1452d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14531795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1454d6ab1dc0SHong Zhang 
1455d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1456d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1457d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1458d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1459d6ab1dc0SHong Zhang   ai    = a_loc->i;
1460d6ab1dc0SHong Zhang   aj    = a_loc->j;
1461d6ab1dc0SHong Zhang 
1462d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1463d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1464d6ab1dc0SHong Zhang     adj = aj + ai[i];
1465d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1466d6ab1dc0SHong Zhang 
1467d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1468e745005bSHong Zhang     /*-------------------------------------------------------------*/
1469d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1470d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1471d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1472d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1473d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1474d6ab1dc0SHong Zhang       row = poJ[j];
1475d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1476d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1477e745005bSHong Zhang       /* perform sparse axpy */
1478e745005bSHong Zhang       nexta  = 0;
1479d6ab1dc0SHong Zhang       valtmp = pA[j];
1480e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1481e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1482e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1483e745005bSHong Zhang           nexta++;
1484d6ab1dc0SHong Zhang         }
1485e745005bSHong Zhang       }
1486e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1487d6ab1dc0SHong Zhang     }
1488d6ab1dc0SHong Zhang 
1489d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1490d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1491d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1492d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1493d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1494d6ab1dc0SHong Zhang       row = pdJ[j];
1495d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1496d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1497e745005bSHong Zhang       /* perform sparse axpy */
1498e745005bSHong Zhang       nexta  = 0;
1499d6ab1dc0SHong Zhang       valtmp = pA[j];
1500e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1501e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1502e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1503e745005bSHong Zhang           nexta++;
1504d6ab1dc0SHong Zhang         }
1505e745005bSHong Zhang       }
1506e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1507d6ab1dc0SHong Zhang     }
1508d6ab1dc0SHong Zhang   }
1509d6ab1dc0SHong Zhang 
1510d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1511d6ab1dc0SHong Zhang   /*------------------------------------*/
1512d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1513d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1514d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1515d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1516d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1517d6ab1dc0SHong Zhang 
1518dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1519d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1520d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1521d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1522d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1523d6ab1dc0SHong Zhang     k++;
1524d6ab1dc0SHong Zhang   }
1525d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1526d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1527d6ab1dc0SHong Zhang 
1528d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1530d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1531d6ab1dc0SHong Zhang 
1532d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1533d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1534dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1535d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1536e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1537d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1538d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1539d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1540d6ab1dc0SHong Zhang   }
1541d6ab1dc0SHong Zhang 
1542d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1543d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1544d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1545d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1546d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1547d6ab1dc0SHong Zhang     /* add received vals into ba */
1548d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1549d6ab1dc0SHong Zhang       /* i-th row */
1550d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1551d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1552d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1553d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1554d6ab1dc0SHong Zhang         nextcj = 0;
1555d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1556d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1557d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1558d6ab1dc0SHong Zhang           }
1559d6ab1dc0SHong Zhang         }
1560d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1561d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1562d6ab1dc0SHong Zhang       }
1563d6ab1dc0SHong Zhang     }
1564d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1565d6ab1dc0SHong Zhang   }
1566d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1567d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1568d6ab1dc0SHong Zhang 
1569d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1570d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1571d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1572d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1573d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1574d6ab1dc0SHong Zhang }
1575d6ab1dc0SHong Zhang 
15762bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
15772bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1578d6ab1dc0SHong Zhang #undef __FUNCT__
15796da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
15806da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1581d6ab1dc0SHong Zhang {
1582d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1583d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1584d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15850298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1586d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1587d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1588d6ab1dc0SHong Zhang   PetscInt            nnz;
1589d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1590d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1591ce94432eSBarry Smith   MPI_Comm            comm;
1592d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1593d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1594d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1595d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1596d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1597d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1598d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1599d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1600d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1601d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1602c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1603d6ab1dc0SHong Zhang   PetscScalar         *vals;
1604d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1605d6ab1dc0SHong Zhang 
1606d6ab1dc0SHong Zhang   PetscFunctionBegin;
1607ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1608d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1609d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1610d6ab1dc0SHong 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);
1611d6ab1dc0SHong Zhang   }
1612d6ab1dc0SHong Zhang 
1613d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1614d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1615d6ab1dc0SHong Zhang 
1616d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1617b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1618d6ab1dc0SHong Zhang 
1619d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1620d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16212205254eSKarl Rupp 
1622d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1623d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1624d6ab1dc0SHong Zhang   ai          = a_loc->i;
1625d6ab1dc0SHong Zhang   aj          = a_loc->j;
1626d6ab1dc0SHong Zhang 
1627d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1628d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1629d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1630d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1631d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1632d6ab1dc0SHong Zhang 
1633d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1634d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1635d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1636d6ab1dc0SHong Zhang 
1637d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1638d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1639d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1640785e854fSJed Brown   ierr   = PetscMalloc1((pon+1),&coi);CHKERRQ(ierr);
1641d6ab1dc0SHong Zhang   coi[0] = 0;
1642d6ab1dc0SHong Zhang 
1643d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1644d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1645a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1646d6ab1dc0SHong Zhang   current_space = free_space;
164719f0e57aSHong Zhang 
1648d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1649d6ab1dc0SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1650c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1651d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1652d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1653d6ab1dc0SHong Zhang 
1654d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1655d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1656d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1657d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1658d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1659d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1660d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1661d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1662d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1663d6ab1dc0SHong Zhang     }
1664d6ab1dc0SHong Zhang     nnz = lnk[0];
1665d6ab1dc0SHong Zhang 
1666d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1667d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1668d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1669d6ab1dc0SHong Zhang       nspacedouble++;
1670d6ab1dc0SHong Zhang     }
1671d6ab1dc0SHong Zhang 
1672d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1673d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16742205254eSKarl Rupp 
1675d6ab1dc0SHong Zhang     current_space->array           += nnz;
1676d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1677d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16782205254eSKarl Rupp 
1679d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1680d6ab1dc0SHong Zhang   }
1681d6ab1dc0SHong Zhang 
1682785e854fSJed Brown   ierr = PetscMalloc1((coi[pon]+1),&coj);CHKERRQ(ierr);
1683d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16842205254eSKarl Rupp 
1685118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1686d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1687d6ab1dc0SHong Zhang 
1688d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1689d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1690d6ab1dc0SHong Zhang   /* determine row ownership */
1691b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1692d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
16932205254eSKarl Rupp 
1694d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1695d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
16962205254eSKarl Rupp 
1697d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1698d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1699d6ab1dc0SHong Zhang 
1700d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17011795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1702785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17032205254eSKarl Rupp 
1704d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1705d6ab1dc0SHong Zhang   merge->nsend = 0;
1706d6ab1dc0SHong Zhang 
1707785e854fSJed Brown   ierr = PetscMalloc1((size+2),&owners_co);CHKERRQ(ierr);
1708d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1709d6ab1dc0SHong Zhang 
1710d6ab1dc0SHong Zhang   proc = 0;
1711d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1712d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1713d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1714d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1715d6ab1dc0SHong Zhang   }
1716d6ab1dc0SHong Zhang 
1717d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1718d6ab1dc0SHong Zhang   owners_co[0] = 0;
1719d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1720d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1721d6ab1dc0SHong Zhang     if (len_si[proc]) {
1722d6ab1dc0SHong Zhang       merge->nsend++;
1723d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1724d6ab1dc0SHong Zhang       len         += len_si[proc];
1725d6ab1dc0SHong Zhang     }
1726d6ab1dc0SHong Zhang   }
1727d6ab1dc0SHong Zhang 
1728d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17290298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1730d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1731d6ab1dc0SHong Zhang 
1732d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1733d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1734d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1735785e854fSJed Brown   ierr = PetscMalloc1((merge->nsend+1),&swaits);CHKERRQ(ierr);
1736d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1737d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1738d6ab1dc0SHong Zhang     i    = owners_co[proc];
1739d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1740d6ab1dc0SHong Zhang     k++;
1741d6ab1dc0SHong Zhang   }
1742d6ab1dc0SHong Zhang 
1743d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1744785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1745d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1746c280ed6aSJed Brown     PetscMPIInt icompleted;
1747c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1748d6ab1dc0SHong Zhang   }
1749d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1750d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1751d6ab1dc0SHong Zhang 
1752d6ab1dc0SHong Zhang   /* send and recv coi */
1753d6ab1dc0SHong Zhang   /*-------------------*/
1754d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1755d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1756785e854fSJed Brown   ierr   = PetscMalloc1((len+1),&buf_s);CHKERRQ(ierr);
1757d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1758d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1759d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1760d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1761d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1762d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1763d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1764d6ab1dc0SHong Zhang     */
1765d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1766d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1767d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1768d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1769d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1770d6ab1dc0SHong Zhang     nrows       = 0;
1771d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1772d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1773d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1774d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1775d6ab1dc0SHong Zhang       nrows++;
1776d6ab1dc0SHong Zhang     }
1777d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1778d6ab1dc0SHong Zhang     k++;
1779d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1780d6ab1dc0SHong Zhang   }
1781d6ab1dc0SHong Zhang   i = merge->nrecv;
1782d6ab1dc0SHong Zhang   while (i--) {
1783c280ed6aSJed Brown     PetscMPIInt icompleted;
1784c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1785d6ab1dc0SHong Zhang   }
1786d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1787d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1788d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1789d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1790d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1791d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1792d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1793d6ab1dc0SHong Zhang 
1794d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1795d6ab1dc0SHong Zhang   /*------------------------------------------*/
1796d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1797785e854fSJed Brown   ierr  = PetscMalloc1((pn+1),&bi);CHKERRQ(ierr);
1798d6ab1dc0SHong Zhang   bi[0] = 0;
1799d6ab1dc0SHong Zhang 
1800d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1801d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1802a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1803d6ab1dc0SHong Zhang   current_space = free_space;
1804d6ab1dc0SHong Zhang 
1805dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1806d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1807d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1808d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1809d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18102205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1811d6ab1dc0SHong Zhang   }
1812d6ab1dc0SHong Zhang 
1813d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1814d6ab1dc0SHong Zhang   rmax = 0;
1815d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1816d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1817d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1818d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1819d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1820d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1821d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1822d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1823d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1824d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1825d6ab1dc0SHong Zhang     }
1826d6ab1dc0SHong Zhang 
1827d6ab1dc0SHong Zhang     /* add received col data into lnk */
1828d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1829d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1830d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1831d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1832d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1833d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1834d6ab1dc0SHong Zhang       }
1835d6ab1dc0SHong Zhang     }
1836d6ab1dc0SHong Zhang     nnz = lnk[0];
1837d6ab1dc0SHong Zhang 
1838d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1839d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1840d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1841d6ab1dc0SHong Zhang       nspacedouble++;
1842d6ab1dc0SHong Zhang     }
1843d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1844d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1845d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18462205254eSKarl Rupp 
1847d6ab1dc0SHong Zhang     current_space->array           += nnz;
1848d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1849d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18502205254eSKarl Rupp 
1851d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1852d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1853d6ab1dc0SHong Zhang   }
1854d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1855d6ab1dc0SHong Zhang 
1856785e854fSJed Brown   ierr      = PetscMalloc1((bi[pn]+1),&bj);CHKERRQ(ierr);
1857d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1858118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1859d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1860d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1861d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1862d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1863d6ab1dc0SHong Zhang 
1864d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1865d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
18661795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1867d6ab1dc0SHong Zhang 
1868d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1869d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
187033d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1871d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1872d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1873d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1874d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1875d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1876d6ab1dc0SHong Zhang     row  = i + rstart;
1877d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1878d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1879d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1880d6ab1dc0SHong Zhang   }
1881d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1882d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1883d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1884d6ab1dc0SHong Zhang 
1885d6ab1dc0SHong Zhang   merge->bi        = bi;
1886d6ab1dc0SHong Zhang   merge->bj        = bj;
1887d6ab1dc0SHong Zhang   merge->coi       = coi;
1888d6ab1dc0SHong Zhang   merge->coj       = coj;
1889d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1890d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1891d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1892d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1893d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1894d6ab1dc0SHong Zhang 
18956da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1896d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1897d6ab1dc0SHong Zhang 
1898d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1899d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19002205254eSKarl Rupp 
1901d6ab1dc0SHong Zhang   c->ptap     = ptap;
19020298fd71SBarry Smith   ptap->api   = NULL;
19030298fd71SBarry Smith   ptap->apj   = NULL;
1904d6ab1dc0SHong Zhang   ptap->merge = merge;
1905d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
19060298fd71SBarry Smith   ptap->apa   = NULL;
1907f96d37fcSHong Zhang 
1908f96d37fcSHong Zhang   *C = Cmpi;
1909f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1910f96d37fcSHong Zhang   if (bi[pn] != 0) {
191157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
191257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1913f96d37fcSHong Zhang   } else {
1914f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1915f96d37fcSHong Zhang   }
1916f96d37fcSHong Zhang #endif
1917f96d37fcSHong Zhang   PetscFunctionReturn(0);
1918f96d37fcSHong Zhang }
1919