xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 2bbb1c24870617abef987903b26f7b4cf1ab8c9e)
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;
100598bc09dSHong Zhang 
101598bc09dSHong Zhang   PetscFunctionBegin;
102a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
103598bc09dSHong Zhang   /*-----------------------------------------------------*/
104a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
105b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
106a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
107598bc09dSHong Zhang 
108598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
109598bc09dSHong Zhang   /*----------------------------------------------------------*/
110598bc09dSHong Zhang   /* get data from symbolic products */
111a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
112a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
113598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
114598bc09dSHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
115598bc09dSHong Zhang 
116598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
117598bc09dSHong Zhang   apa = ptap->apa;
118598bc09dSHong Zhang 
11929825010SHong Zhang   api = ptap->api;
12029825010SHong Zhang   apj = ptap->apj;
121598bc09dSHong Zhang   for (i=0; i<cm; i++) {
122598bc09dSHong Zhang     /* diagonal portion of A */
123598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
124598bc09dSHong Zhang     adj = ad->j + adi[i];
125598bc09dSHong Zhang     ada = ad->a + adi[i];
126598bc09dSHong Zhang     for (j=0; j<anz; j++) {
127598bc09dSHong Zhang       row = adj[j];
128598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
129598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
130598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
131598bc09dSHong Zhang 
132598bc09dSHong Zhang       /* perform dense axpy */
133598bc09dSHong Zhang       valtmp = ada[j];
134598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
135598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
136598bc09dSHong Zhang       }
137598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
138598bc09dSHong Zhang     }
139598bc09dSHong Zhang 
140598bc09dSHong Zhang     /* off-diagonal portion of A */
141598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
142598bc09dSHong Zhang     aoj = ao->j + aoi[i];
143598bc09dSHong Zhang     aoa = ao->a + aoi[i];
144598bc09dSHong Zhang     for (j=0; j<anz; j++) {
145598bc09dSHong Zhang       row = aoj[j];
146598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
147598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
148598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
149598bc09dSHong Zhang 
150598bc09dSHong Zhang       /* perform dense axpy */
151598bc09dSHong Zhang       valtmp = aoa[j];
152598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
153598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
154598bc09dSHong Zhang       }
155598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
156598bc09dSHong Zhang     }
157598bc09dSHong Zhang 
158598bc09dSHong Zhang     /* set values in C */
159598bc09dSHong Zhang     apJ  = apj + api[i];
160598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
161598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
162598bc09dSHong Zhang 
163598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
164598bc09dSHong Zhang     ca = coa + co->i[i];
165598bc09dSHong Zhang     k  = 0;
166598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
167598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
168598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
169598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
170598bc09dSHong Zhang       k++;
171598bc09dSHong Zhang     }
172598bc09dSHong Zhang 
173598bc09dSHong Zhang     /* diagonal part of C */
174598bc09dSHong Zhang     ca = cda + cd->i[i];
175598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
176598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
177598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
178598bc09dSHong Zhang       k++;
179598bc09dSHong Zhang     }
180598bc09dSHong Zhang 
181598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
182598bc09dSHong Zhang     ca = coa + co->i[i];
183598bc09dSHong Zhang     for (; k0<conz; k0++) {
184598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
185598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
186598bc09dSHong Zhang       k++;
187598bc09dSHong Zhang     }
188598bc09dSHong Zhang   }
189598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
190598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191598bc09dSHong Zhang   PetscFunctionReturn(0);
192598bc09dSHong Zhang }
193598bc09dSHong Zhang 
194598bc09dSHong Zhang #undef __FUNCT__
1950fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
1960fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
197598bc09dSHong Zhang {
198598bc09dSHong Zhang   PetscErrorCode     ierr;
199ce94432eSBarry Smith   MPI_Comm           comm;
200bfcd1a73SHong Zhang   Mat                Cmpi;
201598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
2020298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2034ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
204bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2054ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2064ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
207d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
208bfcd1a73SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
209598bc09dSHong Zhang   PetscBT            lnkbt;
210598bc09dSHong Zhang   PetscScalar        *apa;
211bfcd1a73SHong Zhang   PetscReal          afill;
212f84c536bSHong Zhang   PetscInt           nlnk_max,armax,prmax;
213598bc09dSHong Zhang 
214598bc09dSHong Zhang   PetscFunctionBegin;
215ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
216a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
217cf1a0672SHong 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);
218a1a4e84aSHong Zhang   }
219152983bfSHong Zhang 
220598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
221598bc09dSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
222598bc09dSHong Zhang 
223598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
224b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
22519f0e57aSHong Zhang 
226598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
227a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
228598bc09dSHong Zhang 
229a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
230a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
231598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
232598bc09dSHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
233598bc09dSHong Zhang 
234598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
235598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
236598bc09dSHong Zhang   ierr      = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
237a1a4e84aSHong Zhang   ptap->api = api;
238598bc09dSHong Zhang   api[0]    = 0;
239598bc09dSHong Zhang 
240598bc09dSHong Zhang   /* create and initialize a linked list */
241f84c536bSHong Zhang   armax    = ad->rmax+ao->rmax;
242f84c536bSHong Zhang   prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
243f84c536bSHong Zhang   nlnk_max = armax*prmax;
244f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2450d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
246598bc09dSHong Zhang 
247bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
248bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
2492205254eSKarl Rupp 
250598bc09dSHong Zhang   current_space = free_space;
251598bc09dSHong Zhang 
252598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
253598bc09dSHong Zhang   for (i=0; i<am; i++) {
254598bc09dSHong Zhang     /* diagonal portion of A */
255598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
256598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
257598bc09dSHong Zhang       row  = *adj++;
258598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
259598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
260598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2611fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
262598bc09dSHong Zhang     }
263598bc09dSHong Zhang     /* off-diagonal portion of A */
264598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
265598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
266598bc09dSHong Zhang       row  = *aoj++;
267598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
268598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2691fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
270598bc09dSHong Zhang     }
271598bc09dSHong Zhang 
272d14fa243SHong Zhang     apnz     = lnk[0];
273598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
274598bc09dSHong Zhang 
275598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
276598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
277598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
278598bc09dSHong Zhang       nspacedouble++;
279598bc09dSHong Zhang     }
280598bc09dSHong Zhang 
281598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2821fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
283598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2842205254eSKarl Rupp 
285598bc09dSHong Zhang     current_space->array           += apnz;
286598bc09dSHong Zhang     current_space->local_used      += apnz;
287598bc09dSHong Zhang     current_space->local_remaining -= apnz;
288598bc09dSHong Zhang   }
289598bc09dSHong Zhang 
290598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
291598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
292a1a4e84aSHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
293a1a4e84aSHong Zhang   apj  = ptap->apj;
294a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
295598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
296598bc09dSHong Zhang 
297f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
298f84c536bSHong Zhang   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
299f84c536bSHong Zhang   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
3002205254eSKarl Rupp 
301f84c536bSHong Zhang   ptap->apa = apa;
302f84c536bSHong Zhang 
303bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
304598bc09dSHong Zhang   /*----------------------------------------------------*/
305bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
306bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
307a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,A->rmap->bs,P->cmap->bs);CHKERRQ(ierr);
308a2f3521dSMark F. Adams 
309bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
310bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
311598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
312598bc09dSHong Zhang   for (i=0; i<am; i++) {
313598bc09dSHong Zhang     row  = i + rstart;
314598bc09dSHong Zhang     apnz = api[i+1] - api[i];
315bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
316598bc09dSHong Zhang     apj += apnz;
317598bc09dSHong Zhang   }
318bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
319bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
320598bc09dSHong Zhang 
321bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
322bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
323bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
324bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
325598bc09dSHong Zhang 
326bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
327bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
328598bc09dSHong Zhang   c->ptap = ptap;
329598bc09dSHong Zhang 
330bfcd1a73SHong Zhang   *C = Cmpi;
331bfcd1a73SHong Zhang 
332bfcd1a73SHong Zhang   /* set MatInfo */
333118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
334bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
335bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
336bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
337bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
338bfcd1a73SHong Zhang 
339bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
340bfcd1a73SHong Zhang   if (api[am]) {
341bfcd1a73SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
342bfcd1a73SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
343bfcd1a73SHong Zhang   } else {
344bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
345bfcd1a73SHong Zhang   }
346bfcd1a73SHong Zhang #endif
347598bc09dSHong Zhang   PetscFunctionReturn(0);
348598bc09dSHong Zhang }
349598bc09dSHong Zhang 
3509123193aSHong Zhang #undef __FUNCT__
3519123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3529123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3539123193aSHong Zhang {
3549123193aSHong Zhang   PetscErrorCode ierr;
3559123193aSHong Zhang 
3569123193aSHong Zhang   PetscFunctionBegin;
3579123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3583ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3599123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3603ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3619123193aSHong Zhang   }
3623ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3639123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3643ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3659123193aSHong Zhang   PetscFunctionReturn(0);
3669123193aSHong Zhang }
3679123193aSHong Zhang 
3684324174eSBarry Smith typedef struct {
3694324174eSBarry Smith   Mat         workB;
3704324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3714324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3724324174eSBarry Smith } MPIAIJ_MPIDense;
3734324174eSBarry Smith 
3747af9e4fdSHong Zhang #undef __FUNCT__
37596b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
37696b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3774324174eSBarry Smith {
3784324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3794324174eSBarry Smith   PetscErrorCode  ierr;
3804324174eSBarry Smith 
3814324174eSBarry Smith   PetscFunctionBegin;
3826bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3834324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3844324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3854324174eSBarry Smith   PetscFunctionReturn(0);
3864324174eSBarry Smith }
3874324174eSBarry Smith 
3889123193aSHong Zhang #undef __FUNCT__
3899123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
3909123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
3919123193aSHong Zhang {
3929123193aSHong Zhang   PetscErrorCode         ierr;
393f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
394d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
395bf0cc555SLisandro Dalcin   PetscContainer         container;
3964324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
3974324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
3984324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
3994324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
400d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4019123193aSHong Zhang 
4029123193aSHong Zhang   PetscFunctionBegin;
403ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
404d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
405a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(*C,A->rmap->bs,B->cmap->bs);CHKERRQ(ierr);
406cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4070298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
408cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4102205254eSKarl Rupp 
411d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
412f32d5d43SBarry Smith 
4134324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
414f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4150298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4164324174eSBarry Smith   /* Create work arrays needed */
417d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
418d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
4194324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
4204324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
4214324174eSBarry Smith 
422ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
423bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
42496b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
425bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
426bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4279123193aSHong Zhang   PetscFunctionReturn(0);
4289123193aSHong Zhang }
4299123193aSHong Zhang 
4307af9e4fdSHong Zhang #undef __FUNCT__
4317af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
432f32d5d43SBarry Smith /*
4332636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4342636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
435f32d5d43SBarry Smith */
4364324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
437f32d5d43SBarry Smith {
438f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
439f32d5d43SBarry Smith   PetscErrorCode         ierr;
440f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
441f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
442f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
443f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
444f32d5d43SBarry Smith   PetscInt               i,j,k;
445aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
446aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
447f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
448ce94432eSBarry Smith   MPI_Comm               comm;
449d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
450f32d5d43SBarry Smith   MPI_Status             status;
4514324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
452bf0cc555SLisandro Dalcin   PetscContainer         container;
4534324174eSBarry Smith   Mat                    workB;
454f32d5d43SBarry Smith 
455f32d5d43SBarry Smith   PetscFunctionBegin;
456ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
457bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
45829825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
459bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4604324174eSBarry Smith 
4614324174eSBarry Smith   workB = *outworkB = contents->workB;
462cf1a0672SHong 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);
463f32d5d43SBarry Smith   sindices = to->indices;
464f32d5d43SBarry Smith   sstarts  = to->starts;
465f32d5d43SBarry Smith   sprocs   = to->procs;
4664324174eSBarry Smith   swaits   = contents->swaits;
4674324174eSBarry Smith   svalues  = contents->svalues;
468f32d5d43SBarry Smith 
469f32d5d43SBarry Smith   rindices = from->indices;
470f32d5d43SBarry Smith   rstarts  = from->starts;
471f32d5d43SBarry Smith   rprocs   = from->procs;
4724324174eSBarry Smith   rwaits   = contents->rwaits;
4734324174eSBarry Smith   rvalues  = contents->rvalues;
474f32d5d43SBarry Smith 
4758c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
4768c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
477f32d5d43SBarry Smith 
478f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
479f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
480f32d5d43SBarry Smith   }
481f32d5d43SBarry Smith 
482f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
483f32d5d43SBarry Smith     /* pack a message at a time */
484f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
485f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
486f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
4872636ff26SBarry Smith       }
4882636ff26SBarry Smith     }
489f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
490f32d5d43SBarry Smith   }
4912636ff26SBarry Smith 
492f32d5d43SBarry Smith   nrecvs = from->n;
493f32d5d43SBarry Smith   while (nrecvs) {
494f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
495f32d5d43SBarry Smith     nrecvs--;
496f32d5d43SBarry Smith     /* unpack a message at a time */
497f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
498f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
499f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5002636ff26SBarry Smith       }
5012636ff26SBarry Smith     }
502f32d5d43SBarry Smith   }
503cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
504f32d5d43SBarry Smith 
5058c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5068c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
507f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
508f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
509f32d5d43SBarry Smith   PetscFunctionReturn(0);
510f32d5d43SBarry Smith }
511f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
512f32d5d43SBarry Smith 
5139123193aSHong Zhang #undef __FUNCT__
5149123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5159123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5169123193aSHong Zhang {
5179123193aSHong Zhang   PetscErrorCode ierr;
518f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
519f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
520f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
521f32d5d43SBarry Smith   Mat            workB;
5229123193aSHong Zhang 
5239123193aSHong Zhang   PetscFunctionBegin;
524f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
525f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
526f32d5d43SBarry Smith 
527f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5284324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
529f32d5d43SBarry Smith 
530f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
531f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5329123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5339123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5349123193aSHong Zhang   PetscFunctionReturn(0);
5359123193aSHong Zhang }
536cf1a0672SHong Zhang 
5371634b032SHong Zhang #undef __FUNCT__
538cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
539cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
5401634b032SHong Zhang {
5411634b032SHong Zhang   PetscErrorCode ierr;
542cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
543cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
544cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
545cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
546cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
547cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
548cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
54929825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
55029825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
551cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
55229825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
553cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
55429825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
55529825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
5561634b032SHong Zhang 
5571634b032SHong Zhang   PetscFunctionBegin;
558cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
559cf1a0672SHong Zhang   /*-----------------------------------------------------*/
560cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
561b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
562cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
563cf1a0672SHong Zhang 
564cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
565cf1a0672SHong Zhang   /*----------------------------------------------------------*/
566cf1a0672SHong Zhang   /* get data from symbolic products */
567cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
568cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
569cf1a0672SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
570cf1a0672SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
571cf1a0672SHong Zhang 
57229825010SHong Zhang   api = ptap->api;
57329825010SHong Zhang   apj = ptap->apj;
574cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
57529825010SHong Zhang     apJ = apj + api[i];
57629825010SHong Zhang 
577cf1a0672SHong Zhang     /* diagonal portion of A */
578cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
579cf1a0672SHong Zhang     adj = ad->j + adi[i];
580cf1a0672SHong Zhang     ada = ad->a + adi[i];
581cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
582cf1a0672SHong Zhang       row = adj[j];
583cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
584cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
585cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
58629825010SHong Zhang       /* perform sparse axpy */
587cf1a0672SHong Zhang       valtmp = ada[j];
58829825010SHong Zhang       nextp  = 0;
58929825010SHong Zhang       for (k=0; nextp<pnz; k++) {
59029825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
59129825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
59229825010SHong Zhang         }
5931634b032SHong Zhang       }
594cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
595cf1a0672SHong Zhang     }
5961634b032SHong Zhang 
597cf1a0672SHong Zhang     /* off-diagonal portion of A */
598cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
599cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
600cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
601cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
602cf1a0672SHong Zhang       row = aoj[j];
603cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
604cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
605cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
60629825010SHong Zhang       /* perform sparse axpy */
607cf1a0672SHong Zhang       valtmp = aoa[j];
60829825010SHong Zhang       nextp  = 0;
60929825010SHong Zhang       for (k=0; nextp<pnz; k++) {
61029825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
61129825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
61229825010SHong Zhang         }
613cf1a0672SHong Zhang       }
614cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
615cf1a0672SHong Zhang     }
6161634b032SHong Zhang 
617cf1a0672SHong Zhang     /* set values in C */
618cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
619cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6201634b032SHong Zhang 
621cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
622cf1a0672SHong Zhang     ca = coa + co->i[i];
623cf1a0672SHong Zhang     k  = 0;
624cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
625cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
62629825010SHong Zhang       ca[k0]        = apa_sparse[k];
62729825010SHong Zhang       apa_sparse[k] = 0.0;
628cf1a0672SHong Zhang       k++;
629cf1a0672SHong Zhang     }
6301634b032SHong Zhang 
631cf1a0672SHong Zhang     /* diagonal part of C */
632cf1a0672SHong Zhang     ca = cda + cd->i[i];
633cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
63429825010SHong Zhang       ca[k1]        = apa_sparse[k];
63529825010SHong Zhang       apa_sparse[k] = 0.0;
636cf1a0672SHong Zhang       k++;
637cf1a0672SHong Zhang     }
638cf1a0672SHong Zhang 
639cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
640cf1a0672SHong Zhang     ca = coa + co->i[i];
641cf1a0672SHong Zhang     for (; k0<conz; k0++) {
64229825010SHong Zhang       ca[k0]        = apa_sparse[k];
64329825010SHong Zhang       apa_sparse[k] = 0.0;
644cf1a0672SHong Zhang       k++;
645cf1a0672SHong Zhang     }
646cf1a0672SHong Zhang   }
647cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
648cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
649cf1a0672SHong Zhang   PetscFunctionReturn(0);
650cf1a0672SHong Zhang }
651cf1a0672SHong Zhang 
6520fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
653cf1a0672SHong Zhang #undef __FUNCT__
654b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
655b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
656cf1a0672SHong Zhang {
657cf1a0672SHong Zhang   PetscErrorCode     ierr;
658ce94432eSBarry Smith   MPI_Comm           comm;
659cf1a0672SHong Zhang   Mat                Cmpi;
660cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6610298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
662cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
663cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
664cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
665cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
666f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
667cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
668cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
669cf1a0672SHong Zhang   PetscReal          afill;
670cf1a0672SHong Zhang   PetscScalar        *apa;
671cf1a0672SHong Zhang 
672cf1a0672SHong Zhang   PetscFunctionBegin;
673ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
674cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
675cf1a0672SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
676cf1a0672SHong Zhang 
677cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
678b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
67919f0e57aSHong Zhang 
680cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
681cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
682cf1a0672SHong Zhang 
683cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
684cf1a0672SHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
685cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
686cf1a0672SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
687cf1a0672SHong Zhang 
688cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
689cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
690cf1a0672SHong Zhang   ierr      = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
691cf1a0672SHong Zhang   ptap->api = api;
692cf1a0672SHong Zhang   api[0]    = 0;
693cf1a0672SHong Zhang 
694cf1a0672SHong Zhang   /* create and initialize a linked list */
695cf1a0672SHong Zhang   armax    = ad->rmax+ao->rmax;
696cf1a0672SHong Zhang   prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
697cf1a0672SHong Zhang   nlnk_max = armax*prmax;
698cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
699f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
700cf1a0672SHong Zhang 
701cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
702cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
7032205254eSKarl Rupp 
704cf1a0672SHong Zhang   current_space = free_space;
705cf1a0672SHong Zhang 
706cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
707cf1a0672SHong Zhang   for (i=0; i<am; i++) {
708cf1a0672SHong Zhang     /* diagonal portion of A */
709cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
710cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
711cf1a0672SHong Zhang       row  = *adj++;
712cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
713cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
714cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
715f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
716cf1a0672SHong Zhang     }
717cf1a0672SHong Zhang     /* off-diagonal portion of A */
718cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
719cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
720cf1a0672SHong Zhang       row  = *aoj++;
721cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
722cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
723f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
724cf1a0672SHong Zhang     }
725cf1a0672SHong Zhang 
726f84c536bSHong Zhang     apnz     = *lnk;
727cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
728cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
729cf1a0672SHong Zhang 
730cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
731cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
732cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
733cf1a0672SHong Zhang       nspacedouble++;
734cf1a0672SHong Zhang     }
735cf1a0672SHong Zhang 
736cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
737f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
738cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7392205254eSKarl Rupp 
740cf1a0672SHong Zhang     current_space->array           += apnz;
741cf1a0672SHong Zhang     current_space->local_used      += apnz;
742cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
743cf1a0672SHong Zhang   }
744cf1a0672SHong Zhang 
745cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
746cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
747cf1a0672SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
748cf1a0672SHong Zhang   apj  = ptap->apj;
749cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
750f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
751cf1a0672SHong Zhang 
752cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
753cf1a0672SHong Zhang   /*----------------------------------------------------*/
754cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
755cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
756a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,A->rmap->bs,P->cmap->bs);CHKERRQ(ierr);
757cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
758cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
759cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
760cf1a0672SHong Zhang 
76129825010SHong Zhang   /* malloc apa for assembly Cmpi */
762cf1a0672SHong Zhang   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
763cf1a0672SHong Zhang   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
7642205254eSKarl Rupp 
76529825010SHong Zhang   ptap->apa = apa;
766cf1a0672SHong Zhang   for (i=0; i<am; i++) {
767cf1a0672SHong Zhang     row  = i + rstart;
768cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
769cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
770cf1a0672SHong Zhang     apj += apnz;
771cf1a0672SHong Zhang   }
772cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
773cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
774cf1a0672SHong Zhang 
775cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
776cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
777cf1a0672SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
778cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
779cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
780cf1a0672SHong Zhang 
781cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
782cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
783cf1a0672SHong Zhang   c->ptap = ptap;
784cf1a0672SHong Zhang 
785cf1a0672SHong Zhang   *C = Cmpi;
786cf1a0672SHong Zhang 
787cf1a0672SHong Zhang   /* set MatInfo */
788118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
789cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
790cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
791cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
792cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
793cf1a0672SHong Zhang 
794cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
795cf1a0672SHong Zhang   if (api[am]) {
796cf1a0672SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
797cf1a0672SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
798cf1a0672SHong Zhang   } else {
799cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
800cf1a0672SHong Zhang   }
801cf1a0672SHong Zhang #endif
8021634b032SHong Zhang   PetscFunctionReturn(0);
8031634b032SHong Zhang }
804f96d37fcSHong Zhang 
805f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
806f96d37fcSHong Zhang #undef __FUNCT__
807f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
808c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
809f96d37fcSHong Zhang {
810f96d37fcSHong Zhang   PetscErrorCode ierr;
811c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
812c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
813f96d37fcSHong Zhang 
814f96d37fcSHong Zhang   PetscFunctionBegin;
815c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
816c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
817c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
818c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
819c216dbf3SHong Zhang 
820c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
821c216dbf3SHong Zhang     switch (alg) {
822c216dbf3SHong Zhang     case 1:
823*2bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
824c216dbf3SHong Zhang       break;
825c216dbf3SHong Zhang     case 2:
826ab1b013aSHong Zhang       Mat         Pt;
827ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
828ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
829ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
830ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
831ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
832ab1b013aSHong Zhang       ptap     = c->ptap;
833ab1b013aSHong Zhang       ptap->Pt = Pt;
834c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
835c216dbf3SHong Zhang       PetscFunctionReturn(0);
836c216dbf3SHong Zhang       break;
837c216dbf3SHong Zhang     default:
838c216dbf3SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr);
839c216dbf3SHong Zhang       break;
840c216dbf3SHong Zhang     }
841c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
842c216dbf3SHong Zhang   }
843c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
844c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
845c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
846ab1b013aSHong Zhang   PetscFunctionReturn(0);
847ab1b013aSHong Zhang }
848ab1b013aSHong Zhang 
849c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
850c216dbf3SHong Zhang #undef __FUNCT__
851c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
852c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
853c216dbf3SHong Zhang {
854c216dbf3SHong Zhang   PetscErrorCode ierr;
855*2bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
856*2bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
857*2bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
858c216dbf3SHong Zhang 
859c216dbf3SHong Zhang   PetscFunctionBegin;
860c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
861c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
862f96d37fcSHong Zhang   PetscFunctionReturn(0);
863f96d37fcSHong Zhang }
864f96d37fcSHong Zhang 
865f96d37fcSHong Zhang #undef __FUNCT__
866f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
867c5af039cSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
868f96d37fcSHong Zhang {
869c5af039cSHong Zhang   PetscErrorCode      ierr;
870c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
871497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
872c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
873c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
874d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
875497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
876e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
877c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
878ce94432eSBarry Smith   MPI_Comm            comm;
879c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
880c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
881c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
882c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
883c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
884c5af039cSHong Zhang   MPI_Status          *status;
885c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
886d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
887a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
888c5af039cSHong Zhang   Mat                 A_loc;
889c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
890f96d37fcSHong Zhang 
891f96d37fcSHong Zhang   PetscFunctionBegin;
892ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
893c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
894c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
895c5af039cSHong Zhang 
896c5af039cSHong Zhang   ptap  = c->ptap;
897c5af039cSHong Zhang   merge = ptap->merge;
898c5af039cSHong Zhang 
899c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
900c5af039cSHong Zhang   /*--------------------------------------------------------------*/
901c5af039cSHong Zhang   /* get data from symbolic products */
902c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
903c5af039cSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
904c5af039cSHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
905c5af039cSHong Zhang 
906c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
907c5af039cSHong Zhang   owners = merge->rowmap->range;
908c5af039cSHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
909c5af039cSHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
910c5af039cSHong Zhang 
911c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
912c5af039cSHong Zhang   A_loc = ptap->A_loc;
913c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
914c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
915d6ab1dc0SHong Zhang   ai    = a_loc->i;
916d6ab1dc0SHong Zhang   aj    = a_loc->j;
917c5af039cSHong Zhang 
918e745005bSHong Zhang   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */
919e745005bSHong Zhang   ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
920c5af039cSHong Zhang 
921c5af039cSHong Zhang   for (i=0; i<am; i++) {
922e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
923d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
924d6ab1dc0SHong Zhang     adj = aj + ai[i];
925d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
926c5af039cSHong Zhang     for (j=0; j<anz; j++) {
927e745005bSHong Zhang       aval[adj[j]] = ada[j];
928c5af039cSHong Zhang     }
929c5af039cSHong Zhang 
930c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
931c5af039cSHong Zhang     /*--------------------------------------------------------------*/
932c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
933c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
934c5af039cSHong Zhang     poJ = po->j + po->i[i];
935c5af039cSHong Zhang     pA  = po->a + po->i[i];
936c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
937c5af039cSHong Zhang       row = poJ[j];
938c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
939c5af039cSHong Zhang       cj  = coj + coi[row];
940c5af039cSHong Zhang       ca  = coa + coi[row];
941c5af039cSHong Zhang       /* perform dense axpy */
942c5af039cSHong Zhang       valtmp = pA[j];
943c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
944e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
945c5af039cSHong Zhang       }
946c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
947c5af039cSHong Zhang     }
948c5af039cSHong Zhang 
949c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
950c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
951c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
952c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
953c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
954c5af039cSHong Zhang       row = pdJ[j];
955c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
956c5af039cSHong Zhang       cj  = bj + bi[row];
957c5af039cSHong Zhang       ca  = ba + bi[row];
958c5af039cSHong Zhang       /* perform dense axpy */
959c5af039cSHong Zhang       valtmp = pA[j];
960c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
961e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
962c5af039cSHong Zhang       }
963c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
964c5af039cSHong Zhang     }
965c5af039cSHong Zhang 
966d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
967d6ab1dc0SHong Zhang     aJ = aj + ai[i];
968e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
969c5af039cSHong Zhang   }
970c5af039cSHong Zhang 
971c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
972c5af039cSHong Zhang   /*------------------------------------*/
973c5af039cSHong Zhang   buf_ri = merge->buf_ri;
974c5af039cSHong Zhang   buf_rj = merge->buf_rj;
975c5af039cSHong Zhang   len_s  = merge->len_s;
976c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
977c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
978c5af039cSHong Zhang 
979c5af039cSHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
980c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
981c5af039cSHong Zhang     if (!len_s[proc]) continue;
982c5af039cSHong Zhang     i    = merge->owners_co[proc];
983c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
984c5af039cSHong Zhang     k++;
985c5af039cSHong Zhang   }
986c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
987c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
988c5af039cSHong Zhang 
989c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
990c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
991c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
992c5af039cSHong Zhang 
993c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
994c5af039cSHong Zhang   /*----------------------------------------------------*/
995c5af039cSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
996c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
997c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
998c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
999c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1000c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1001c5af039cSHong Zhang   }
1002c5af039cSHong Zhang 
1003c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1004c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1005c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1006c5af039cSHong Zhang     ba_i = ba + bi[i];
1007c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1008c5af039cSHong Zhang     /* add received vals into ba */
1009c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1010c5af039cSHong Zhang       /* i-th row */
1011c5af039cSHong Zhang       if (i == *nextrow[k]) {
1012c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1013c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1014c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1015c5af039cSHong Zhang         nextcj = 0;
1016c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1017c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1018c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1019c5af039cSHong Zhang           }
1020c5af039cSHong Zhang         }
1021c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1022c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1023c5af039cSHong Zhang       }
1024c5af039cSHong Zhang     }
1025c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1026c5af039cSHong Zhang   }
1027c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1028c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1029c5af039cSHong Zhang 
1030c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1031c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1032c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1033c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1034e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1035f96d37fcSHong Zhang   PetscFunctionReturn(0);
1036f96d37fcSHong Zhang }
1037f96d37fcSHong Zhang 
1038c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1039f96d37fcSHong Zhang #undef __FUNCT__
1040*2bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
1041*2bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1042f96d37fcSHong Zhang {
1043f96d37fcSHong Zhang   PetscErrorCode      ierr;
10444e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1045f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10460298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1047f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1048f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1049f96d37fcSHong Zhang   PetscInt            nnz;
10504e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1051497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1052f96d37fcSHong Zhang   PetscBT             lnkbt;
1053ce94432eSBarry Smith   MPI_Comm            comm;
1054f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1055f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1056f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1057f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1058f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1059f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1060f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1061f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1062d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1063f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1064c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1065f96d37fcSHong Zhang   PetscScalar         *vals;
10664e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1067f96d37fcSHong Zhang 
1068f96d37fcSHong Zhang   PetscFunctionBegin;
1069ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1070c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1071c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1072c5af039cSHong 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);
1073c5af039cSHong Zhang   }
1074c5af039cSHong Zhang 
1075f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1076f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1077f96d37fcSHong Zhang 
1078f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1079f96d37fcSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1080f96d37fcSHong Zhang 
1081f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1082f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
10832205254eSKarl Rupp 
1084c5af039cSHong Zhang   ptap->A_loc = A_loc;
10852205254eSKarl Rupp 
10861c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1087d6ab1dc0SHong Zhang   ai    = a_loc->i;
1088d6ab1dc0SHong Zhang   aj    = a_loc->j;
1089f96d37fcSHong Zhang 
1090f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1091f96d37fcSHong Zhang   /*----------------------------------------------------*/
10924e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
10934e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
10944e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
10954e938277SHong Zhang 
10964e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
10974e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
10984e938277SHong Zhang   poti = pot->i; potj = pot->j;
1099f96d37fcSHong Zhang 
1100f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11012205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1102f96d37fcSHong Zhang   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1103f96d37fcSHong Zhang   coi[0] = 0;
1104f96d37fcSHong Zhang 
1105f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1106d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1107a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1108f96d37fcSHong Zhang   current_space = free_space;
110919f0e57aSHong Zhang 
1110f96d37fcSHong Zhang   /* create and initialize a linked list */
11114e938277SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1112c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
11134e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
11144e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1115f96d37fcSHong Zhang 
1116f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1117f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1118f96d37fcSHong Zhang     ptJ = potj + poti[i];
1119f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1120f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1121d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1122d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1123f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1124d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1125f96d37fcSHong Zhang     }
11264e938277SHong Zhang     nnz = lnk[0];
1127f96d37fcSHong Zhang 
1128f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1129f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1130f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1131f96d37fcSHong Zhang       nspacedouble++;
1132f96d37fcSHong Zhang     }
1133f96d37fcSHong Zhang 
1134f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11354e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11362205254eSKarl Rupp 
1137f96d37fcSHong Zhang     current_space->array           += nnz;
1138f96d37fcSHong Zhang     current_space->local_used      += nnz;
1139f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11402205254eSKarl Rupp 
1141f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1142f96d37fcSHong Zhang   }
1143f96d37fcSHong Zhang 
1144f96d37fcSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1145f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11462205254eSKarl Rupp 
1147118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1148f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1149f96d37fcSHong Zhang 
1150f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1151f96d37fcSHong Zhang   /*----------------------------------------------*/
1152f96d37fcSHong Zhang   /* determine row ownership */
1153f96d37fcSHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1154f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11552205254eSKarl Rupp 
1156f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1157f96d37fcSHong Zhang   merge->rowmap->bs = 1;
11582205254eSKarl Rupp 
1159f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1160f96d37fcSHong Zhang   owners = merge->rowmap->range;
1161f96d37fcSHong Zhang 
1162f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
1163f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1164f96d37fcSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1165f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
11662205254eSKarl Rupp 
1167f96d37fcSHong Zhang   len_s        = merge->len_s;
1168f96d37fcSHong Zhang   merge->nsend = 0;
1169f96d37fcSHong Zhang 
1170f96d37fcSHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1171f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1172f96d37fcSHong Zhang 
1173f96d37fcSHong Zhang   proc = 0;
1174f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1175f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1176f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1177f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1178f96d37fcSHong Zhang   }
1179f96d37fcSHong Zhang 
1180f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1181f96d37fcSHong Zhang   owners_co[0] = 0;
1182f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1183f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1184f96d37fcSHong Zhang     if (len_si[proc]) {
1185f96d37fcSHong Zhang       merge->nsend++;
1186f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1187f96d37fcSHong Zhang       len         += len_si[proc];
1188f96d37fcSHong Zhang     }
1189f96d37fcSHong Zhang   }
1190f96d37fcSHong Zhang 
1191f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
11920298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1193f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1194f96d37fcSHong Zhang 
1195f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1196f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1197f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1198f96d37fcSHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1199f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1200f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1201f96d37fcSHong Zhang     i    = owners_co[proc];
1202f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1203f96d37fcSHong Zhang     k++;
1204f96d37fcSHong Zhang   }
1205f96d37fcSHong Zhang 
1206f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1207f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1208f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1209c280ed6aSJed Brown     PetscMPIInt icompleted;
1210c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1211f96d37fcSHong Zhang   }
1212f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1213f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1214f96d37fcSHong Zhang 
1215f96d37fcSHong Zhang   /* send and recv coi */
1216f96d37fcSHong Zhang   /*-------------------*/
1217f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1218f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1219f96d37fcSHong Zhang   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1220f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1221f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1222f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1223f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1224f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1225f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1226f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1227f96d37fcSHong Zhang     */
1228f96d37fcSHong Zhang     /*-------------------------------------------*/
1229f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1230f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1231f96d37fcSHong Zhang     buf_si[0]   = nrows;
1232f96d37fcSHong Zhang     buf_si_i[0] = 0;
1233f96d37fcSHong Zhang     nrows       = 0;
1234f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1235f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1236f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1237f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1238f96d37fcSHong Zhang       nrows++;
1239f96d37fcSHong Zhang     }
1240f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1241f96d37fcSHong Zhang     k++;
1242f96d37fcSHong Zhang     buf_si += len_si[proc];
1243f96d37fcSHong Zhang   }
1244f96d37fcSHong Zhang   i = merge->nrecv;
1245f96d37fcSHong Zhang   while (i--) {
1246c280ed6aSJed Brown     PetscMPIInt icompleted;
1247c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1248f96d37fcSHong Zhang   }
1249f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1250f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1251f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1252f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1253f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1254f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1255f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1256f96d37fcSHong Zhang 
1257f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1258f96d37fcSHong Zhang   /*------------------------------------------*/
1259f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1260f96d37fcSHong Zhang   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1261f96d37fcSHong Zhang   bi[0] = 0;
1262f96d37fcSHong Zhang 
1263c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1264d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1265a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1266f96d37fcSHong Zhang   current_space = free_space;
1267f96d37fcSHong Zhang 
1268f96d37fcSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1269f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1270f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1271f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1272f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1273f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1274f96d37fcSHong Zhang   }
1275f96d37fcSHong Zhang 
12761c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1277f96d37fcSHong Zhang   rmax = 0;
1278f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1279f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1280f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1281f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1282f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1283f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1284d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1285d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1286f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1287d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1288f96d37fcSHong Zhang     }
12894e938277SHong Zhang 
1290f96d37fcSHong Zhang     /* add received col data into lnk */
1291f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1292f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1293f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1294f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
12954e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1296f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1297f96d37fcSHong Zhang       }
1298f96d37fcSHong Zhang     }
12994e938277SHong Zhang     nnz = lnk[0];
1300f96d37fcSHong Zhang 
1301f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1302f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1303f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1304f96d37fcSHong Zhang       nspacedouble++;
1305f96d37fcSHong Zhang     }
1306f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13074e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1308f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13092205254eSKarl Rupp 
1310f96d37fcSHong Zhang     current_space->array           += nnz;
1311f96d37fcSHong Zhang     current_space->local_used      += nnz;
1312f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13132205254eSKarl Rupp 
1314f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1315f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1316f96d37fcSHong Zhang   }
1317f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1318f96d37fcSHong Zhang 
1319f96d37fcSHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1320f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13212205254eSKarl Rupp 
1322118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1323f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1324d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13254e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13264e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1327f96d37fcSHong Zhang 
13281c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13291c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13301c7d5954SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
13311c7d5954SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1332f96d37fcSHong Zhang 
1333f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1334f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1335a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1336f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1337f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1338f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1339f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1340f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1341f96d37fcSHong Zhang     row  = i + rstart;
1342f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1343f96d37fcSHong Zhang     Jptr = bj + bi[i];
1344f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1345f96d37fcSHong Zhang   }
1346f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1347f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1348f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1349f96d37fcSHong Zhang 
1350f96d37fcSHong Zhang   merge->bi        = bi;
1351f96d37fcSHong Zhang   merge->bj        = bj;
1352f96d37fcSHong Zhang   merge->coi       = coi;
1353f96d37fcSHong Zhang   merge->coj       = coj;
1354f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1355f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1356f96d37fcSHong Zhang   merge->owners_co = owners_co;
1357f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1358f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1359f96d37fcSHong Zhang 
1360d6ab1dc0SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1361c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1362f96d37fcSHong Zhang 
1363f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1364f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1365f96d37fcSHong Zhang   c->ptap     = ptap;
13660298fd71SBarry Smith   ptap->api   = NULL;
13670298fd71SBarry Smith   ptap->apj   = NULL;
1368f96d37fcSHong Zhang   ptap->merge = merge;
1369d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1370d6ab1dc0SHong Zhang 
1371d6ab1dc0SHong Zhang   *C = Cmpi;
1372d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1373d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
1374d6ab1dc0SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1375d6ab1dc0SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1376d6ab1dc0SHong Zhang   } else {
1377d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1378d6ab1dc0SHong Zhang   }
1379d6ab1dc0SHong Zhang #endif
1380d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1381d6ab1dc0SHong Zhang }
1382d6ab1dc0SHong Zhang 
1383d6ab1dc0SHong Zhang #undef __FUNCT__
1384d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
1385d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C)
1386d6ab1dc0SHong Zhang {
1387d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1388d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1389d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1390d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1391d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1392e745005bSHong Zhang   PetscInt            *adj;
1393e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1394e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1395d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1396ce94432eSBarry Smith   MPI_Comm            comm;
1397d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1398d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1399d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1400d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1401d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1402d6ab1dc0SHong Zhang   MPI_Status          *status;
1403d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1404d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1405a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1406d6ab1dc0SHong Zhang   Mat                 A_loc;
1407d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1408d6ab1dc0SHong Zhang 
1409d6ab1dc0SHong Zhang   PetscFunctionBegin;
1410ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1411d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1412d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1413d6ab1dc0SHong Zhang 
1414d6ab1dc0SHong Zhang   ptap  = c->ptap;
1415d6ab1dc0SHong Zhang   merge = ptap->merge;
1416d6ab1dc0SHong Zhang 
1417e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1418e745005bSHong Zhang   /*------------------------------------------*/
1419d6ab1dc0SHong Zhang   /* get data from symbolic products */
1420d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1421d6ab1dc0SHong Zhang   ierr   = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1422d6ab1dc0SHong Zhang   ierr   = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1423d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1424d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1425d6ab1dc0SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1426d6ab1dc0SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1427d6ab1dc0SHong Zhang 
1428d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1429d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1430d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1431d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1432d6ab1dc0SHong Zhang   ai    = a_loc->i;
1433d6ab1dc0SHong Zhang   aj    = a_loc->j;
1434d6ab1dc0SHong Zhang 
1435d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1436d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1437d6ab1dc0SHong Zhang     adj = aj + ai[i];
1438d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1439d6ab1dc0SHong Zhang 
1440d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1441e745005bSHong Zhang     /*-------------------------------------------------------------*/
1442d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1443d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1444d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1445d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1446d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1447d6ab1dc0SHong Zhang       row = poJ[j];
1448d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1449d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1450e745005bSHong Zhang       /* perform sparse axpy */
1451e745005bSHong Zhang       nexta  = 0;
1452d6ab1dc0SHong Zhang       valtmp = pA[j];
1453e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1454e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1455e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1456e745005bSHong Zhang           nexta++;
1457d6ab1dc0SHong Zhang         }
1458e745005bSHong Zhang       }
1459e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1460d6ab1dc0SHong Zhang     }
1461d6ab1dc0SHong Zhang 
1462d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1463d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1464d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1465d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1466d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1467d6ab1dc0SHong Zhang       row = pdJ[j];
1468d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1469d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1470e745005bSHong Zhang       /* perform sparse axpy */
1471e745005bSHong Zhang       nexta  = 0;
1472d6ab1dc0SHong Zhang       valtmp = pA[j];
1473e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1474e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1475e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1476e745005bSHong Zhang           nexta++;
1477d6ab1dc0SHong Zhang         }
1478e745005bSHong Zhang       }
1479e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1480d6ab1dc0SHong Zhang     }
1481d6ab1dc0SHong Zhang   }
1482d6ab1dc0SHong Zhang 
1483d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1484d6ab1dc0SHong Zhang   /*------------------------------------*/
1485d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1486d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1487d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1488d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1489d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1490d6ab1dc0SHong Zhang 
1491d6ab1dc0SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1492d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1493d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1494d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1495d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1496d6ab1dc0SHong Zhang     k++;
1497d6ab1dc0SHong Zhang   }
1498d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1499d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1500d6ab1dc0SHong Zhang 
1501d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1502d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1503d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1504d6ab1dc0SHong Zhang 
1505d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1506d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1507d6ab1dc0SHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1508d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1509e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1510d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1511d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1512d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1513d6ab1dc0SHong Zhang   }
1514d6ab1dc0SHong Zhang 
1515d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1516d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1517d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1518d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1519d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1520d6ab1dc0SHong Zhang     /* add received vals into ba */
1521d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1522d6ab1dc0SHong Zhang       /* i-th row */
1523d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1524d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1525d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1526d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1527d6ab1dc0SHong Zhang         nextcj = 0;
1528d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1529d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1530d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1531d6ab1dc0SHong Zhang           }
1532d6ab1dc0SHong Zhang         }
1533d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1534d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1535d6ab1dc0SHong Zhang       }
1536d6ab1dc0SHong Zhang     }
1537d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1538d6ab1dc0SHong Zhang   }
1539d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1540d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1541d6ab1dc0SHong Zhang 
1542d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1543d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1544d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1545d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1546d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1547d6ab1dc0SHong Zhang }
1548d6ab1dc0SHong Zhang 
1549*2bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
1550*2bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1551d6ab1dc0SHong Zhang #undef __FUNCT__
1552d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
1553d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C)
1554d6ab1dc0SHong Zhang {
1555d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1556d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1557d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15580298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1559d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1560d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1561d6ab1dc0SHong Zhang   PetscInt            nnz;
1562d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1563d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1564ce94432eSBarry Smith   MPI_Comm            comm;
1565d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1566d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1567d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1568d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1569d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1570d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1571d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1572d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1573d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1574d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1575c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1576d6ab1dc0SHong Zhang   PetscScalar         *vals;
1577d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1578d6ab1dc0SHong Zhang 
1579d6ab1dc0SHong Zhang   PetscFunctionBegin;
1580ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1581d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1582d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1583d6ab1dc0SHong 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);
1584d6ab1dc0SHong Zhang   }
1585d6ab1dc0SHong Zhang 
1586d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1587d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1588d6ab1dc0SHong Zhang 
1589d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1590d6ab1dc0SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1591d6ab1dc0SHong Zhang 
1592d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1593d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
15942205254eSKarl Rupp 
1595d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1596d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1597d6ab1dc0SHong Zhang   ai          = a_loc->i;
1598d6ab1dc0SHong Zhang   aj          = a_loc->j;
1599d6ab1dc0SHong Zhang 
1600d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1601d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1602d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1603d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1604d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1605d6ab1dc0SHong Zhang 
1606d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1607d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1608d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1609d6ab1dc0SHong Zhang 
1610d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1611d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1612d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1613d6ab1dc0SHong Zhang   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1614d6ab1dc0SHong Zhang   coi[0] = 0;
1615d6ab1dc0SHong Zhang 
1616d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1617d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1618a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1619d6ab1dc0SHong Zhang   current_space = free_space;
162019f0e57aSHong Zhang 
1621d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1622d6ab1dc0SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1623c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1624d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1625d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1626d6ab1dc0SHong Zhang 
1627d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1628d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1629d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1630d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1631d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1632d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1633d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1634d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1635d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1636d6ab1dc0SHong Zhang     }
1637d6ab1dc0SHong Zhang     nnz = lnk[0];
1638d6ab1dc0SHong Zhang 
1639d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1640d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1641d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1642d6ab1dc0SHong Zhang       nspacedouble++;
1643d6ab1dc0SHong Zhang     }
1644d6ab1dc0SHong Zhang 
1645d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1646d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16472205254eSKarl Rupp 
1648d6ab1dc0SHong Zhang     current_space->array           += nnz;
1649d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1650d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16512205254eSKarl Rupp 
1652d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1653d6ab1dc0SHong Zhang   }
1654d6ab1dc0SHong Zhang 
1655d6ab1dc0SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1656d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16572205254eSKarl Rupp 
1658118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1659d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1660d6ab1dc0SHong Zhang 
1661d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1662d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1663d6ab1dc0SHong Zhang   /* determine row ownership */
1664d6ab1dc0SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1665d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
16662205254eSKarl Rupp 
1667d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1668d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
16692205254eSKarl Rupp 
1670d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1671d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1672d6ab1dc0SHong Zhang 
1673d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
1674d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1675d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1676d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
16772205254eSKarl Rupp 
1678d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1679d6ab1dc0SHong Zhang   merge->nsend = 0;
1680d6ab1dc0SHong Zhang 
1681d6ab1dc0SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1682d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1683d6ab1dc0SHong Zhang 
1684d6ab1dc0SHong Zhang   proc = 0;
1685d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1686d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1687d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1688d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1689d6ab1dc0SHong Zhang   }
1690d6ab1dc0SHong Zhang 
1691d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1692d6ab1dc0SHong Zhang   owners_co[0] = 0;
1693d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1694d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1695d6ab1dc0SHong Zhang     if (len_si[proc]) {
1696d6ab1dc0SHong Zhang       merge->nsend++;
1697d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1698d6ab1dc0SHong Zhang       len         += len_si[proc];
1699d6ab1dc0SHong Zhang     }
1700d6ab1dc0SHong Zhang   }
1701d6ab1dc0SHong Zhang 
1702d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17030298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1704d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1705d6ab1dc0SHong Zhang 
1706d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1707d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1708d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1709d6ab1dc0SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1710d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1711d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1712d6ab1dc0SHong Zhang     i    = owners_co[proc];
1713d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1714d6ab1dc0SHong Zhang     k++;
1715d6ab1dc0SHong Zhang   }
1716d6ab1dc0SHong Zhang 
1717d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1718d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1719d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1720c280ed6aSJed Brown     PetscMPIInt icompleted;
1721c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1722d6ab1dc0SHong Zhang   }
1723d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1724d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1725d6ab1dc0SHong Zhang 
1726d6ab1dc0SHong Zhang   /* send and recv coi */
1727d6ab1dc0SHong Zhang   /*-------------------*/
1728d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1729d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1730d6ab1dc0SHong Zhang   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1731d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1732d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1733d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1734d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1735d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1736d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1737d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1738d6ab1dc0SHong Zhang     */
1739d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1740d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1741d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1742d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1743d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1744d6ab1dc0SHong Zhang     nrows       = 0;
1745d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1746d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1747d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1748d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1749d6ab1dc0SHong Zhang       nrows++;
1750d6ab1dc0SHong Zhang     }
1751d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1752d6ab1dc0SHong Zhang     k++;
1753d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1754d6ab1dc0SHong Zhang   }
1755d6ab1dc0SHong Zhang   i = merge->nrecv;
1756d6ab1dc0SHong Zhang   while (i--) {
1757c280ed6aSJed Brown     PetscMPIInt icompleted;
1758c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1759d6ab1dc0SHong Zhang   }
1760d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1761d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1762d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1763d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1764d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1765d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1766d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1767d6ab1dc0SHong Zhang 
1768d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1769d6ab1dc0SHong Zhang   /*------------------------------------------*/
1770d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1771d6ab1dc0SHong Zhang   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1772d6ab1dc0SHong Zhang   bi[0] = 0;
1773d6ab1dc0SHong Zhang 
1774d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1775d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1776a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1777d6ab1dc0SHong Zhang   current_space = free_space;
1778d6ab1dc0SHong Zhang 
1779d6ab1dc0SHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1780d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1781d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1782d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1783d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
17842205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1785d6ab1dc0SHong Zhang   }
1786d6ab1dc0SHong Zhang 
1787d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1788d6ab1dc0SHong Zhang   rmax = 0;
1789d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1790d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1791d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1792d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1793d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1794d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1795d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1796d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1797d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1798d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1799d6ab1dc0SHong Zhang     }
1800d6ab1dc0SHong Zhang 
1801d6ab1dc0SHong Zhang     /* add received col data into lnk */
1802d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1803d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1804d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1805d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1806d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1807d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1808d6ab1dc0SHong Zhang       }
1809d6ab1dc0SHong Zhang     }
1810d6ab1dc0SHong Zhang     nnz = lnk[0];
1811d6ab1dc0SHong Zhang 
1812d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1813d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1814d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1815d6ab1dc0SHong Zhang       nspacedouble++;
1816d6ab1dc0SHong Zhang     }
1817d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1818d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18202205254eSKarl Rupp 
1821d6ab1dc0SHong Zhang     current_space->array           += nnz;
1822d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1823d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18242205254eSKarl Rupp 
1825d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1826d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1827d6ab1dc0SHong Zhang   }
1828d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1829d6ab1dc0SHong Zhang 
1830d6ab1dc0SHong Zhang   ierr      = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1831d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1832118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1833d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1834d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1835d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1836d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1837d6ab1dc0SHong Zhang 
1838d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1839d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
1840d6ab1dc0SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1841d6ab1dc0SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1842d6ab1dc0SHong Zhang 
1843d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1844d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1845a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1846d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1847d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1848d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1849d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1850d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1851d6ab1dc0SHong Zhang     row  = i + rstart;
1852d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1853d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1854d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1855d6ab1dc0SHong Zhang   }
1856d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1857d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1858d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1859d6ab1dc0SHong Zhang 
1860d6ab1dc0SHong Zhang   merge->bi        = bi;
1861d6ab1dc0SHong Zhang   merge->bj        = bj;
1862d6ab1dc0SHong Zhang   merge->coi       = coi;
1863d6ab1dc0SHong Zhang   merge->coj       = coj;
1864d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1865d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1866d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1867d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1868d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1869d6ab1dc0SHong Zhang 
1870d6ab1dc0SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1871d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1872d6ab1dc0SHong Zhang 
1873d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1874d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
18752205254eSKarl Rupp 
1876d6ab1dc0SHong Zhang   c->ptap     = ptap;
18770298fd71SBarry Smith   ptap->api   = NULL;
18780298fd71SBarry Smith   ptap->apj   = NULL;
1879d6ab1dc0SHong Zhang   ptap->merge = merge;
1880d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
18810298fd71SBarry Smith   ptap->apa   = NULL;
1882f96d37fcSHong Zhang 
1883f96d37fcSHong Zhang   *C = Cmpi;
1884f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1885f96d37fcSHong Zhang   if (bi[pn] != 0) {
1886f96d37fcSHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
18871c7d5954SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1888f96d37fcSHong Zhang   } else {
1889f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1890f96d37fcSHong Zhang   }
1891f96d37fcSHong Zhang #endif
1892f96d37fcSHong Zhang   PetscFunctionReturn(0);
1893f96d37fcSHong Zhang }
1894