xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision c9ba551fbcb486bcce22a477c8127a5b0fe113f5)
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__
83f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
84f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(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;
1009467f45fSHong Zhang   MPI_Comm       comm;
1019467f45fSHong Zhang   PetscMPIInt    size;
102598bc09dSHong Zhang 
103598bc09dSHong Zhang   PetscFunctionBegin;
1049467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1059467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1069467f45fSHong Zhang 
107a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
108598bc09dSHong Zhang   /*-----------------------------------------------------*/
109a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
110b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
111a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
112598bc09dSHong Zhang 
113598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
114598bc09dSHong Zhang   /*----------------------------------------------------------*/
115598bc09dSHong Zhang   /* get data from symbolic products */
116a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
117598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
1189467f45fSHong Zhang   if (size >1) {
1199467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
120598bc09dSHong Zhang     pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
121968382fdSHong Zhang   } else {
122968382fdSHong Zhang     pi_oth=NULL; pj_oth=NULL; pa_oth=NULL;
1239467f45fSHong Zhang   }
124598bc09dSHong Zhang 
125598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
126598bc09dSHong Zhang   apa = ptap->apa;
127598bc09dSHong Zhang 
12829825010SHong Zhang   api = ptap->api;
12929825010SHong Zhang   apj = ptap->apj;
130598bc09dSHong Zhang   for (i=0; i<cm; i++) {
131598bc09dSHong Zhang     /* diagonal portion of A */
132598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
133598bc09dSHong Zhang     adj = ad->j + adi[i];
134598bc09dSHong Zhang     ada = ad->a + adi[i];
135598bc09dSHong Zhang     for (j=0; j<anz; j++) {
136598bc09dSHong Zhang       row = adj[j];
137598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
138598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
139598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
140598bc09dSHong Zhang 
141598bc09dSHong Zhang       /* perform dense axpy */
142598bc09dSHong Zhang       valtmp = ada[j];
143598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
144598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
145598bc09dSHong Zhang       }
146598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
147598bc09dSHong Zhang     }
148598bc09dSHong Zhang 
149598bc09dSHong Zhang     /* off-diagonal portion of A */
150598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
151598bc09dSHong Zhang     aoj = ao->j + aoi[i];
152598bc09dSHong Zhang     aoa = ao->a + aoi[i];
153598bc09dSHong Zhang     for (j=0; j<anz; j++) {
154598bc09dSHong Zhang       row = aoj[j];
155598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
156598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
157598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
158598bc09dSHong Zhang 
159598bc09dSHong Zhang       /* perform dense axpy */
160598bc09dSHong Zhang       valtmp = aoa[j];
161598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
162598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
163598bc09dSHong Zhang       }
164598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
165598bc09dSHong Zhang     }
166598bc09dSHong Zhang 
167598bc09dSHong Zhang     /* set values in C */
168598bc09dSHong Zhang     apJ  = apj + api[i];
169598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
170598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
171598bc09dSHong Zhang 
172598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
173598bc09dSHong Zhang     ca = coa + co->i[i];
174598bc09dSHong Zhang     k  = 0;
175598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
176598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
177598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
178598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
179598bc09dSHong Zhang       k++;
180598bc09dSHong Zhang     }
181598bc09dSHong Zhang 
182598bc09dSHong Zhang     /* diagonal part of C */
183598bc09dSHong Zhang     ca = cda + cd->i[i];
184598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
185598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
186598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
187598bc09dSHong Zhang       k++;
188598bc09dSHong Zhang     }
189598bc09dSHong Zhang 
190598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
191598bc09dSHong Zhang     ca = coa + co->i[i];
192598bc09dSHong Zhang     for (; k0<conz; k0++) {
193598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
194598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
195598bc09dSHong Zhang       k++;
196598bc09dSHong Zhang     }
197598bc09dSHong Zhang   }
198598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
200598bc09dSHong Zhang   PetscFunctionReturn(0);
201598bc09dSHong Zhang }
202598bc09dSHong Zhang 
203598bc09dSHong Zhang #undef __FUNCT__
2040fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
2050fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
206598bc09dSHong Zhang {
207598bc09dSHong Zhang   PetscErrorCode     ierr;
208ce94432eSBarry Smith   MPI_Comm           comm;
2099467f45fSHong Zhang   PetscMPIInt        size;
210bfcd1a73SHong Zhang   Mat                Cmpi;
211598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
2120298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2134ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
214bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2154ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2164ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
217d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
218bfcd1a73SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
219598bc09dSHong Zhang   PetscBT            lnkbt;
220598bc09dSHong Zhang   PetscScalar        *apa;
221bfcd1a73SHong Zhang   PetscReal          afill;
222f84c536bSHong Zhang   PetscInt           nlnk_max,armax,prmax;
223598bc09dSHong Zhang 
224598bc09dSHong Zhang   PetscFunctionBegin;
225ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
2269467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2279467f45fSHong Zhang 
228a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
229cf1a0672SHong 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);
230a1a4e84aSHong Zhang   }
231152983bfSHong Zhang 
232598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
233b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
234598bc09dSHong Zhang 
235598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
236b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
23719f0e57aSHong Zhang 
238598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
239a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
240598bc09dSHong Zhang 
241a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
242598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2439467f45fSHong Zhang   if (size > 1) {
2449467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
245598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
24620e3dc75SHong Zhang   } else {
24720e3dc75SHong Zhang     p_oth = NULL;
24820e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2499467f45fSHong Zhang   }
250598bc09dSHong Zhang 
251598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
252598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
253854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
254a1a4e84aSHong Zhang   ptap->api = api;
255598bc09dSHong Zhang   api[0]    = 0;
256598bc09dSHong Zhang 
257598bc09dSHong Zhang   /* create and initialize a linked list */
258f84c536bSHong Zhang   armax    = ad->rmax+ao->rmax;
2599467f45fSHong Zhang   if (size >1) {
260f84c536bSHong Zhang     prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
2619467f45fSHong Zhang   } else {
2629467f45fSHong Zhang     prmax = p_loc->rmax;
2639467f45fSHong Zhang   }
264f84c536bSHong Zhang   nlnk_max = armax*prmax;
265f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2660d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
267598bc09dSHong Zhang 
268bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
269bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
2702205254eSKarl Rupp 
271598bc09dSHong Zhang   current_space = free_space;
272598bc09dSHong Zhang 
273598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
274598bc09dSHong Zhang   for (i=0; i<am; i++) {
275598bc09dSHong Zhang     /* diagonal portion of A */
276598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
277598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
278598bc09dSHong Zhang       row  = *adj++;
279598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
280598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
281598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2821fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
283598bc09dSHong Zhang     }
284598bc09dSHong Zhang     /* off-diagonal portion of A */
285598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
286598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
287598bc09dSHong Zhang       row  = *aoj++;
288598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
289598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2901fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
291598bc09dSHong Zhang     }
292598bc09dSHong Zhang 
293d14fa243SHong Zhang     apnz     = lnk[0];
294598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
295598bc09dSHong Zhang 
296598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
297598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
298598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
299598bc09dSHong Zhang       nspacedouble++;
300598bc09dSHong Zhang     }
301598bc09dSHong Zhang 
302598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
3031fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
304598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
3052205254eSKarl Rupp 
306598bc09dSHong Zhang     current_space->array           += apnz;
307598bc09dSHong Zhang     current_space->local_used      += apnz;
308598bc09dSHong Zhang     current_space->local_remaining -= apnz;
309598bc09dSHong Zhang   }
310598bc09dSHong Zhang 
311598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
312598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
313854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
314a1a4e84aSHong Zhang   apj  = ptap->apj;
315a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
316598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
317598bc09dSHong Zhang 
318f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
3191795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
3202205254eSKarl Rupp 
321f84c536bSHong Zhang   ptap->apa = apa;
322f84c536bSHong Zhang 
323bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
324598bc09dSHong Zhang   /*----------------------------------------------------*/
325bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
326bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
32733d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
328a2f3521dSMark F. Adams 
329bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
330bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
331598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
332598bc09dSHong Zhang   for (i=0; i<am; i++) {
333598bc09dSHong Zhang     row  = i + rstart;
334598bc09dSHong Zhang     apnz = api[i+1] - api[i];
335bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
336598bc09dSHong Zhang     apj += apnz;
337598bc09dSHong Zhang   }
338bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
339bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
340598bc09dSHong Zhang 
341bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
342bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
343f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
344bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
345bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
346598bc09dSHong Zhang 
347bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
348bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
349598bc09dSHong Zhang   c->ptap = ptap;
350598bc09dSHong Zhang 
351bfcd1a73SHong Zhang   *C = Cmpi;
352bfcd1a73SHong Zhang 
353bfcd1a73SHong Zhang   /* set MatInfo */
354118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
355bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
356bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
357bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
358bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
359bfcd1a73SHong Zhang 
360bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
361bfcd1a73SHong Zhang   if (api[am]) {
36257622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
36357622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
364bfcd1a73SHong Zhang   } else {
365bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
366bfcd1a73SHong Zhang   }
367bfcd1a73SHong Zhang #endif
368598bc09dSHong Zhang   PetscFunctionReturn(0);
369598bc09dSHong Zhang }
370598bc09dSHong Zhang 
3719123193aSHong Zhang #undef __FUNCT__
3729123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3739123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3749123193aSHong Zhang {
3759123193aSHong Zhang   PetscErrorCode ierr;
3769123193aSHong Zhang 
3779123193aSHong Zhang   PetscFunctionBegin;
3789123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3793ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3809123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3813ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3829123193aSHong Zhang   }
3833ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3849123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3853ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3869123193aSHong Zhang   PetscFunctionReturn(0);
3879123193aSHong Zhang }
3889123193aSHong Zhang 
3894324174eSBarry Smith typedef struct {
3904324174eSBarry Smith   Mat         workB;
3914324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3924324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3934324174eSBarry Smith } MPIAIJ_MPIDense;
3944324174eSBarry Smith 
3957af9e4fdSHong Zhang #undef __FUNCT__
39696b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
39796b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3984324174eSBarry Smith {
3994324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
4004324174eSBarry Smith   PetscErrorCode  ierr;
4014324174eSBarry Smith 
4024324174eSBarry Smith   PetscFunctionBegin;
4036bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
4044324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
4054324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
4064324174eSBarry Smith   PetscFunctionReturn(0);
4074324174eSBarry Smith }
4084324174eSBarry Smith 
4099123193aSHong Zhang #undef __FUNCT__
4109123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4119123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4129123193aSHong Zhang {
4139123193aSHong Zhang   PetscErrorCode         ierr;
414f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
415d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
416bf0cc555SLisandro Dalcin   PetscContainer         container;
4174324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4184324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4194324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4204324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
421d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4229123193aSHong Zhang 
4239123193aSHong Zhang   PetscFunctionBegin;
424ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
425d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
42633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
427cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4280298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
429cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
430cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4312205254eSKarl Rupp 
432d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
433f32d5d43SBarry Smith 
434b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
435f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4360298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4374324174eSBarry Smith   /* Create work arrays needed */
438dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
439dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
440dcca6d9dSJed Brown                       from->n,&contents->rwaits,
441dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4424324174eSBarry Smith 
443ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
444bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
44596b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
446bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
447bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4489123193aSHong Zhang   PetscFunctionReturn(0);
4499123193aSHong Zhang }
4509123193aSHong Zhang 
4517af9e4fdSHong Zhang #undef __FUNCT__
4527af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
453f32d5d43SBarry Smith /*
4542636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4552636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
456f32d5d43SBarry Smith */
4574324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
458f32d5d43SBarry Smith {
459f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
460f32d5d43SBarry Smith   PetscErrorCode         ierr;
461f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
462f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
463f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
464f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
465f32d5d43SBarry Smith   PetscInt               i,j,k;
466aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
467aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
468f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
469ce94432eSBarry Smith   MPI_Comm               comm;
470d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
471f32d5d43SBarry Smith   MPI_Status             status;
4724324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
473bf0cc555SLisandro Dalcin   PetscContainer         container;
4744324174eSBarry Smith   Mat                    workB;
475f32d5d43SBarry Smith 
476f32d5d43SBarry Smith   PetscFunctionBegin;
477ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
478bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
47929825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
480bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4814324174eSBarry Smith 
4824324174eSBarry Smith   workB = *outworkB = contents->workB;
483cf1a0672SHong 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);
484f32d5d43SBarry Smith   sindices = to->indices;
485f32d5d43SBarry Smith   sstarts  = to->starts;
486f32d5d43SBarry Smith   sprocs   = to->procs;
4874324174eSBarry Smith   swaits   = contents->swaits;
4884324174eSBarry Smith   svalues  = contents->svalues;
489f32d5d43SBarry Smith 
490f32d5d43SBarry Smith   rindices = from->indices;
491f32d5d43SBarry Smith   rstarts  = from->starts;
492f32d5d43SBarry Smith   rprocs   = from->procs;
4934324174eSBarry Smith   rwaits   = contents->rwaits;
4944324174eSBarry Smith   rvalues  = contents->rvalues;
495f32d5d43SBarry Smith 
4968c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
4978c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
498f32d5d43SBarry Smith 
499f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
500f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
501f32d5d43SBarry Smith   }
502f32d5d43SBarry Smith 
503f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
504f32d5d43SBarry Smith     /* pack a message at a time */
505f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
506f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
507f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5082636ff26SBarry Smith       }
5092636ff26SBarry Smith     }
510f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
511f32d5d43SBarry Smith   }
5122636ff26SBarry Smith 
513f32d5d43SBarry Smith   nrecvs = from->n;
514f32d5d43SBarry Smith   while (nrecvs) {
515f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
516f32d5d43SBarry Smith     nrecvs--;
517f32d5d43SBarry Smith     /* unpack a message at a time */
518f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
519f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
520f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5212636ff26SBarry Smith       }
5222636ff26SBarry Smith     }
523f32d5d43SBarry Smith   }
524cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
525f32d5d43SBarry Smith 
5268c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5278c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
528f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
529f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
530f32d5d43SBarry Smith   PetscFunctionReturn(0);
531f32d5d43SBarry Smith }
532f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
533f32d5d43SBarry Smith 
5349123193aSHong Zhang #undef __FUNCT__
5359123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5369123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5379123193aSHong Zhang {
5389123193aSHong Zhang   PetscErrorCode ierr;
539f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
540f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
541f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
542f32d5d43SBarry Smith   Mat            workB;
5439123193aSHong Zhang 
5449123193aSHong Zhang   PetscFunctionBegin;
545f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
546f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
547f32d5d43SBarry Smith 
548f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5494324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
550f32d5d43SBarry Smith 
551f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
552f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5539123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5549123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5559123193aSHong Zhang   PetscFunctionReturn(0);
5569123193aSHong Zhang }
557cf1a0672SHong Zhang 
5581634b032SHong Zhang #undef __FUNCT__
559f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
560f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5611634b032SHong Zhang {
5621634b032SHong Zhang   PetscErrorCode ierr;
563cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
564cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
565cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
566cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
567cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
568cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
569cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
57029825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
57129825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
572cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
57329825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
574cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
57529825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
57629825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
577a7c7454dSHong Zhang   MPI_Comm       comm;
578a7c7454dSHong Zhang   PetscMPIInt    size;
5791634b032SHong Zhang 
5801634b032SHong Zhang   PetscFunctionBegin;
581a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
582a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
583a7c7454dSHong Zhang 
584cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
585cf1a0672SHong Zhang   /*-----------------------------------------------------*/
586cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
587b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
588cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
589cf1a0672SHong Zhang 
590cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
591cf1a0672SHong Zhang   /*----------------------------------------------------------*/
592cf1a0672SHong Zhang   /* get data from symbolic products */
593cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
594cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
595a7c7454dSHong Zhang   if (size >1) {
596a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
597cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
59820e3dc75SHong Zhang   } else {
59920e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
600a7c7454dSHong Zhang   }
601cf1a0672SHong Zhang 
60229825010SHong Zhang   api = ptap->api;
60329825010SHong Zhang   apj = ptap->apj;
604cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
60529825010SHong Zhang     apJ = apj + api[i];
60629825010SHong Zhang 
607cf1a0672SHong Zhang     /* diagonal portion of A */
608cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
609cf1a0672SHong Zhang     adj = ad->j + adi[i];
610cf1a0672SHong Zhang     ada = ad->a + adi[i];
611cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
612cf1a0672SHong Zhang       row = adj[j];
613cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
614cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
615cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
61629825010SHong Zhang       /* perform sparse axpy */
617cf1a0672SHong Zhang       valtmp = ada[j];
61829825010SHong Zhang       nextp  = 0;
61929825010SHong Zhang       for (k=0; nextp<pnz; k++) {
62029825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
62129825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
62229825010SHong Zhang         }
6231634b032SHong Zhang       }
624cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
625cf1a0672SHong Zhang     }
6261634b032SHong Zhang 
627cf1a0672SHong Zhang     /* off-diagonal portion of A */
628cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
629cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
630cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
631cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
632cf1a0672SHong Zhang       row = aoj[j];
633cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
634cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
635cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
63629825010SHong Zhang       /* perform sparse axpy */
637cf1a0672SHong Zhang       valtmp = aoa[j];
63829825010SHong Zhang       nextp  = 0;
63929825010SHong Zhang       for (k=0; nextp<pnz; k++) {
64029825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
64129825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
64229825010SHong Zhang         }
643cf1a0672SHong Zhang       }
644cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
645cf1a0672SHong Zhang     }
6461634b032SHong Zhang 
647cf1a0672SHong Zhang     /* set values in C */
648cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
649cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6501634b032SHong Zhang 
651cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
652cf1a0672SHong Zhang     ca = coa + co->i[i];
653cf1a0672SHong Zhang     k  = 0;
654cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
655cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
65629825010SHong Zhang       ca[k0]        = apa_sparse[k];
65729825010SHong Zhang       apa_sparse[k] = 0.0;
658cf1a0672SHong Zhang       k++;
659cf1a0672SHong Zhang     }
6601634b032SHong Zhang 
661cf1a0672SHong Zhang     /* diagonal part of C */
662cf1a0672SHong Zhang     ca = cda + cd->i[i];
663cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
66429825010SHong Zhang       ca[k1]        = apa_sparse[k];
66529825010SHong Zhang       apa_sparse[k] = 0.0;
666cf1a0672SHong Zhang       k++;
667cf1a0672SHong Zhang     }
668cf1a0672SHong Zhang 
669cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
670cf1a0672SHong Zhang     ca = coa + co->i[i];
671cf1a0672SHong Zhang     for (; k0<conz; k0++) {
67229825010SHong Zhang       ca[k0]        = apa_sparse[k];
67329825010SHong Zhang       apa_sparse[k] = 0.0;
674cf1a0672SHong Zhang       k++;
675cf1a0672SHong Zhang     }
676cf1a0672SHong Zhang   }
677cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
678cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
679cf1a0672SHong Zhang   PetscFunctionReturn(0);
680cf1a0672SHong Zhang }
681cf1a0672SHong Zhang 
6820fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
683cf1a0672SHong Zhang #undef __FUNCT__
684b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
685b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
686cf1a0672SHong Zhang {
687cf1a0672SHong Zhang   PetscErrorCode     ierr;
688ce94432eSBarry Smith   MPI_Comm           comm;
689a7c7454dSHong Zhang   PetscMPIInt        size;
690cf1a0672SHong Zhang   Mat                Cmpi;
691cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6920298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
693cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
694cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
695cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
696cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
697f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
698cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
699cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
700cf1a0672SHong Zhang   PetscReal          afill;
701cf1a0672SHong Zhang   PetscScalar        *apa;
702cf1a0672SHong Zhang 
703cf1a0672SHong Zhang   PetscFunctionBegin;
704ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
705a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
706a7c7454dSHong Zhang 
707cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
708b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
709cf1a0672SHong Zhang 
710cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
711b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
71219f0e57aSHong Zhang 
713cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
714cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
715cf1a0672SHong Zhang 
716cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
717cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
718a7c7454dSHong Zhang   if (size > 1) {
719a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
72020e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
721968382fdSHong Zhang   } else {
72220e3dc75SHong Zhang     p_oth  = NULL;
72320e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
724a7c7454dSHong Zhang   }
725cf1a0672SHong Zhang 
726cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
727cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
728854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
729cf1a0672SHong Zhang   ptap->api = api;
730cf1a0672SHong Zhang   api[0]    = 0;
731cf1a0672SHong Zhang 
732cf1a0672SHong Zhang   /* create and initialize a linked list */
733cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
734a7c7454dSHong Zhang   if (size >1) {
735cf1a0672SHong Zhang     prmax = PetscMax(p_loc->rmax,p_oth->rmax);
736a7c7454dSHong Zhang   } else {
737a7c7454dSHong Zhang     prmax = p_loc->rmax;
738a7c7454dSHong Zhang   }
739cf1a0672SHong Zhang   nlnk_max = armax*prmax;
740cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
741f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
742cf1a0672SHong Zhang 
743cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
744cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
7452205254eSKarl Rupp 
746cf1a0672SHong Zhang   current_space = free_space;
747cf1a0672SHong Zhang 
748cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
749cf1a0672SHong Zhang   for (i=0; i<am; i++) {
750cf1a0672SHong Zhang     /* diagonal portion of A */
751cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
752cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
753cf1a0672SHong Zhang       row  = *adj++;
754cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
755cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
756cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
757f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
758cf1a0672SHong Zhang     }
759cf1a0672SHong Zhang     /* off-diagonal portion of A */
760cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
761cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
762cf1a0672SHong Zhang       row  = *aoj++;
763cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
764cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
765f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
766cf1a0672SHong Zhang     }
767cf1a0672SHong Zhang 
768f84c536bSHong Zhang     apnz     = *lnk;
769cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
770cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
771cf1a0672SHong Zhang 
772cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
773cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
774cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
775cf1a0672SHong Zhang       nspacedouble++;
776cf1a0672SHong Zhang     }
777cf1a0672SHong Zhang 
778cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
779f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
780cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7812205254eSKarl Rupp 
782cf1a0672SHong Zhang     current_space->array           += apnz;
783cf1a0672SHong Zhang     current_space->local_used      += apnz;
784cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
785cf1a0672SHong Zhang   }
786cf1a0672SHong Zhang 
787cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
788cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
789854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
790cf1a0672SHong Zhang   apj  = ptap->apj;
791cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
792f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
793cf1a0672SHong Zhang 
794cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
795cf1a0672SHong Zhang   /*----------------------------------------------------*/
796cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
797cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
79833d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
799cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
800cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
801cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
802cf1a0672SHong Zhang 
80329825010SHong Zhang   /* malloc apa for assembly Cmpi */
8041795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8052205254eSKarl Rupp 
80629825010SHong Zhang   ptap->apa = apa;
807cf1a0672SHong Zhang   for (i=0; i<am; i++) {
808cf1a0672SHong Zhang     row  = i + rstart;
809cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
810cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
811cf1a0672SHong Zhang     apj += apnz;
812cf1a0672SHong Zhang   }
813cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
814cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
815cf1a0672SHong Zhang 
816cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
817cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
818f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
819cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
820cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
821cf1a0672SHong Zhang 
822cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
823cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
824cf1a0672SHong Zhang   c->ptap = ptap;
825cf1a0672SHong Zhang 
826cf1a0672SHong Zhang   *C = Cmpi;
827cf1a0672SHong Zhang 
828cf1a0672SHong Zhang   /* set MatInfo */
829118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
830cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
831cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
832cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
833cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
834cf1a0672SHong Zhang 
835cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
836cf1a0672SHong Zhang   if (api[am]) {
83757622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
83857622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
839cf1a0672SHong Zhang   } else {
840cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
841cf1a0672SHong Zhang   }
842cf1a0672SHong Zhang #endif
8431634b032SHong Zhang   PetscFunctionReturn(0);
8441634b032SHong Zhang }
845f96d37fcSHong Zhang 
846f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
847f96d37fcSHong Zhang #undef __FUNCT__
848f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
849c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
850f96d37fcSHong Zhang {
851f96d37fcSHong Zhang   PetscErrorCode ierr;
852c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
853c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
854f96d37fcSHong Zhang 
855f96d37fcSHong Zhang   PetscFunctionBegin;
856c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
857c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
858c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
859c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
860c216dbf3SHong Zhang 
861c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
862c216dbf3SHong Zhang     switch (alg) {
863c216dbf3SHong Zhang     case 1:
8642bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
865c216dbf3SHong Zhang       break;
866c216dbf3SHong Zhang     case 2:
867275476c6SMatthew G. Knepley     {
868ab1b013aSHong Zhang       Mat         Pt;
869ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
870ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
871ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
872ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
873ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
874ab1b013aSHong Zhang       ptap     = c->ptap;
875ab1b013aSHong Zhang       ptap->Pt = Pt;
876c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
877c216dbf3SHong Zhang       PetscFunctionReturn(0);
878275476c6SMatthew G. Knepley     }
879c216dbf3SHong Zhang       break;
880c216dbf3SHong Zhang     default:
8816da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
882c216dbf3SHong Zhang       break;
883c216dbf3SHong Zhang     }
884c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
885c216dbf3SHong Zhang   }
886c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
887c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
888c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
889ab1b013aSHong Zhang   PetscFunctionReturn(0);
890ab1b013aSHong Zhang }
891ab1b013aSHong Zhang 
892c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
893c216dbf3SHong Zhang #undef __FUNCT__
894c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
895c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
896c216dbf3SHong Zhang {
897c216dbf3SHong Zhang   PetscErrorCode ierr;
8982bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
8992bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9002bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
901c216dbf3SHong Zhang 
902c216dbf3SHong Zhang   PetscFunctionBegin;
903c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
904c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
905f96d37fcSHong Zhang   PetscFunctionReturn(0);
906f96d37fcSHong Zhang }
907f96d37fcSHong Zhang 
9086da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
909f96d37fcSHong Zhang #undef __FUNCT__
9106da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9116da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
912f96d37fcSHong Zhang {
913c5af039cSHong Zhang   PetscErrorCode      ierr;
914c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
915497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
916c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
917c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
918d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
919497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
920e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
921c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
922ce94432eSBarry Smith   MPI_Comm            comm;
923c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
924c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
925c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
926c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
927c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
928c5af039cSHong Zhang   MPI_Status          *status;
929c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
930d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
931a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
932c5af039cSHong Zhang   Mat                 A_loc;
933c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
934f96d37fcSHong Zhang 
935f96d37fcSHong Zhang   PetscFunctionBegin;
936ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
937c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
938c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
939c5af039cSHong Zhang 
940c5af039cSHong Zhang   ptap  = c->ptap;
941c5af039cSHong Zhang   merge = ptap->merge;
942c5af039cSHong Zhang 
943c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
944c5af039cSHong Zhang   /*--------------------------------------------------------------*/
945c5af039cSHong Zhang   /* get data from symbolic products */
946c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
947854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
948c5af039cSHong Zhang 
949c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
950c5af039cSHong Zhang   owners = merge->rowmap->range;
951854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
952c5af039cSHong Zhang 
953c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
954c5af039cSHong Zhang   A_loc = ptap->A_loc;
955c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
956c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
957d6ab1dc0SHong Zhang   ai    = a_loc->i;
958d6ab1dc0SHong Zhang   aj    = a_loc->j;
959c5af039cSHong Zhang 
960854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
961c5af039cSHong Zhang 
962c5af039cSHong Zhang   for (i=0; i<am; i++) {
963e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
964d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
965d6ab1dc0SHong Zhang     adj = aj + ai[i];
966d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
967c5af039cSHong Zhang     for (j=0; j<anz; j++) {
968e745005bSHong Zhang       aval[adj[j]] = ada[j];
969c5af039cSHong Zhang     }
970c5af039cSHong Zhang 
971c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
972c5af039cSHong Zhang     /*--------------------------------------------------------------*/
973c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
974c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
975c5af039cSHong Zhang     poJ = po->j + po->i[i];
976c5af039cSHong Zhang     pA  = po->a + po->i[i];
977c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
978c5af039cSHong Zhang       row = poJ[j];
979c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
980c5af039cSHong Zhang       cj  = coj + coi[row];
981c5af039cSHong Zhang       ca  = coa + coi[row];
982c5af039cSHong Zhang       /* perform dense axpy */
983c5af039cSHong Zhang       valtmp = pA[j];
984c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
985e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
986c5af039cSHong Zhang       }
987c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
988c5af039cSHong Zhang     }
989c5af039cSHong Zhang 
990c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
991c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
992c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
993c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
994c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
995c5af039cSHong Zhang       row = pdJ[j];
996c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
997c5af039cSHong Zhang       cj  = bj + bi[row];
998c5af039cSHong Zhang       ca  = ba + bi[row];
999c5af039cSHong Zhang       /* perform dense axpy */
1000c5af039cSHong Zhang       valtmp = pA[j];
1001c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1002e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1003c5af039cSHong Zhang       }
1004c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1005c5af039cSHong Zhang     }
1006c5af039cSHong Zhang 
1007d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1008d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1009e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1010c5af039cSHong Zhang   }
1011c5af039cSHong Zhang 
1012c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1013c5af039cSHong Zhang   /*------------------------------------*/
1014c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1015c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1016c5af039cSHong Zhang   len_s  = merge->len_s;
1017c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1018c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1019c5af039cSHong Zhang 
1020dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1021c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1022c5af039cSHong Zhang     if (!len_s[proc]) continue;
1023c5af039cSHong Zhang     i    = merge->owners_co[proc];
1024c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1025c5af039cSHong Zhang     k++;
1026c5af039cSHong Zhang   }
1027c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1028c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1029c5af039cSHong Zhang 
1030c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1031c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1032c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1033c5af039cSHong Zhang 
1034c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1035c5af039cSHong Zhang   /*----------------------------------------------------*/
1036dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1037c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1038c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1039c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1040c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1041c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1042c5af039cSHong Zhang   }
1043c5af039cSHong Zhang 
1044c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1045c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1046c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1047c5af039cSHong Zhang     ba_i = ba + bi[i];
1048c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1049c5af039cSHong Zhang     /* add received vals into ba */
1050c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1051c5af039cSHong Zhang       /* i-th row */
1052c5af039cSHong Zhang       if (i == *nextrow[k]) {
1053c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1054c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1055c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1056c5af039cSHong Zhang         nextcj = 0;
1057c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1058c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1059c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1060c5af039cSHong Zhang           }
1061c5af039cSHong Zhang         }
1062c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1063c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1064c5af039cSHong Zhang       }
1065c5af039cSHong Zhang     }
1066c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1067c5af039cSHong Zhang   }
1068c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1069c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1070c5af039cSHong Zhang 
1071c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1072c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1073c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1074c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1075e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1076f96d37fcSHong Zhang   PetscFunctionReturn(0);
1077f96d37fcSHong Zhang }
1078f96d37fcSHong Zhang 
1079*c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1080c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1081f96d37fcSHong Zhang #undef __FUNCT__
10822bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
10832bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1084f96d37fcSHong Zhang {
1085f96d37fcSHong Zhang   PetscErrorCode      ierr;
10864e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1087f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10880298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1089f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1090f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1091f96d37fcSHong Zhang   PetscInt            nnz;
10924e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1093497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1094f96d37fcSHong Zhang   PetscBT             lnkbt;
1095ce94432eSBarry Smith   MPI_Comm            comm;
1096f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1097f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1098f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1099f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1100f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1101f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1102f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1103f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1104d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1105f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1106c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1107f96d37fcSHong Zhang   PetscScalar         *vals;
11084e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1109f96d37fcSHong Zhang 
1110f96d37fcSHong Zhang   PetscFunctionBegin;
1111ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1112c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1113c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1114c5af039cSHong 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);
1115c5af039cSHong Zhang   }
1116c5af039cSHong Zhang 
1117f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1118f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1119f96d37fcSHong Zhang 
1120f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1121b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1122f96d37fcSHong Zhang 
1123f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1124f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11252205254eSKarl Rupp 
1126c5af039cSHong Zhang   ptap->A_loc = A_loc;
11272205254eSKarl Rupp 
11281c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1129d6ab1dc0SHong Zhang   ai    = a_loc->i;
1130d6ab1dc0SHong Zhang   aj    = a_loc->j;
1131f96d37fcSHong Zhang 
1132f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1133f96d37fcSHong Zhang   /*----------------------------------------------------*/
11344e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11354e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11364e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11374e938277SHong Zhang 
11384e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11394e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11404e938277SHong Zhang   poti = pot->i; potj = pot->j;
1141f96d37fcSHong Zhang 
1142f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11432205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1144854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1145f96d37fcSHong Zhang   coi[0] = 0;
1146f96d37fcSHong Zhang 
1147f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1148d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1149a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1150f96d37fcSHong Zhang   current_space = free_space;
115119f0e57aSHong Zhang 
1152f96d37fcSHong Zhang   /* create and initialize a linked list */
11534e938277SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1154c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
11554e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
11564e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1157f96d37fcSHong Zhang 
1158f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1159f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1160f96d37fcSHong Zhang     ptJ = potj + poti[i];
1161f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1162f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1163d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1164d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1165f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1166d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1167f96d37fcSHong Zhang     }
11684e938277SHong Zhang     nnz = lnk[0];
1169f96d37fcSHong Zhang 
1170f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1171f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1172f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1173f96d37fcSHong Zhang       nspacedouble++;
1174f96d37fcSHong Zhang     }
1175f96d37fcSHong Zhang 
1176f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11774e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11782205254eSKarl Rupp 
1179f96d37fcSHong Zhang     current_space->array           += nnz;
1180f96d37fcSHong Zhang     current_space->local_used      += nnz;
1181f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11822205254eSKarl Rupp 
1183f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1184f96d37fcSHong Zhang   }
1185f96d37fcSHong Zhang 
1186854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1187f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11882205254eSKarl Rupp 
1189118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1190f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1191f96d37fcSHong Zhang 
1192f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1193f96d37fcSHong Zhang   /*----------------------------------------------*/
1194f96d37fcSHong Zhang   /* determine row ownership */
1195b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1196f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11972205254eSKarl Rupp 
1198f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1199f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12002205254eSKarl Rupp 
1201f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1202f96d37fcSHong Zhang   owners = merge->rowmap->range;
1203f96d37fcSHong Zhang 
1204f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12051795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1206785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12072205254eSKarl Rupp 
1208f96d37fcSHong Zhang   len_s        = merge->len_s;
1209f96d37fcSHong Zhang   merge->nsend = 0;
1210f96d37fcSHong Zhang 
1211854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1212f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1213f96d37fcSHong Zhang 
1214f96d37fcSHong Zhang   proc = 0;
1215f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1216f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1217f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1218f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1219f96d37fcSHong Zhang   }
1220f96d37fcSHong Zhang 
1221f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1222f96d37fcSHong Zhang   owners_co[0] = 0;
1223f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1224f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1225f96d37fcSHong Zhang     if (len_si[proc]) {
1226f96d37fcSHong Zhang       merge->nsend++;
1227f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1228f96d37fcSHong Zhang       len         += len_si[proc];
1229f96d37fcSHong Zhang     }
1230f96d37fcSHong Zhang   }
1231f96d37fcSHong Zhang 
1232f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12330298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1234f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1235f96d37fcSHong Zhang 
1236f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1237f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1238f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1239854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1240f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1241f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1242f96d37fcSHong Zhang     i    = owners_co[proc];
1243f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1244f96d37fcSHong Zhang     k++;
1245f96d37fcSHong Zhang   }
1246f96d37fcSHong Zhang 
1247f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1248785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1249f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1250c280ed6aSJed Brown     PetscMPIInt icompleted;
1251c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1252f96d37fcSHong Zhang   }
1253f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1254f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1255f96d37fcSHong Zhang 
1256f96d37fcSHong Zhang   /* send and recv coi */
1257f96d37fcSHong Zhang   /*-------------------*/
1258f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1259f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1260854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1261f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1262f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1263f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1264f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1265f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1266f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1267f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1268f96d37fcSHong Zhang     */
1269f96d37fcSHong Zhang     /*-------------------------------------------*/
1270f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1271f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1272f96d37fcSHong Zhang     buf_si[0]   = nrows;
1273f96d37fcSHong Zhang     buf_si_i[0] = 0;
1274f96d37fcSHong Zhang     nrows       = 0;
1275f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1276f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1277f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1278f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1279f96d37fcSHong Zhang       nrows++;
1280f96d37fcSHong Zhang     }
1281f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1282f96d37fcSHong Zhang     k++;
1283f96d37fcSHong Zhang     buf_si += len_si[proc];
1284f96d37fcSHong Zhang   }
1285f96d37fcSHong Zhang   i = merge->nrecv;
1286f96d37fcSHong Zhang   while (i--) {
1287c280ed6aSJed Brown     PetscMPIInt icompleted;
1288c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1289f96d37fcSHong Zhang   }
1290f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1291f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1292f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1293f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1294f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1295f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1296f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1297f96d37fcSHong Zhang 
1298f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1299f96d37fcSHong Zhang   /*------------------------------------------*/
1300f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1301854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1302f96d37fcSHong Zhang   bi[0] = 0;
1303f96d37fcSHong Zhang 
1304c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1305d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1306a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1307f96d37fcSHong Zhang   current_space = free_space;
1308f96d37fcSHong Zhang 
1309dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1310f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1311f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1312f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1313f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1314f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1315f96d37fcSHong Zhang   }
1316f96d37fcSHong Zhang 
13171c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1318f96d37fcSHong Zhang   rmax = 0;
1319f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1320f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1321f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1322f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1323f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1324f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1325d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1326d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1327f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1328d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1329f96d37fcSHong Zhang     }
13304e938277SHong Zhang 
1331f96d37fcSHong Zhang     /* add received col data into lnk */
1332f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1333f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1334f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1335f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13364e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1337f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1338f96d37fcSHong Zhang       }
1339f96d37fcSHong Zhang     }
13404e938277SHong Zhang     nnz = lnk[0];
1341f96d37fcSHong Zhang 
1342f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1343f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1344f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1345f96d37fcSHong Zhang       nspacedouble++;
1346f96d37fcSHong Zhang     }
1347f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13484e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1349f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13502205254eSKarl Rupp 
1351f96d37fcSHong Zhang     current_space->array           += nnz;
1352f96d37fcSHong Zhang     current_space->local_used      += nnz;
1353f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13542205254eSKarl Rupp 
1355f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1356f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1357f96d37fcSHong Zhang   }
1358f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1359f96d37fcSHong Zhang 
1360854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1361f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13622205254eSKarl Rupp 
1363118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1364f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1365d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13664e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13674e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1368f96d37fcSHong Zhang 
13691c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13701c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13711795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1372f96d37fcSHong Zhang 
1373f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1374f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
137533d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1376f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1377f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1378f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1379f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1380f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1381f96d37fcSHong Zhang     row  = i + rstart;
1382f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1383f96d37fcSHong Zhang     Jptr = bj + bi[i];
1384f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1385f96d37fcSHong Zhang   }
1386f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1387f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1388f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1389f96d37fcSHong Zhang 
1390f96d37fcSHong Zhang   merge->bi        = bi;
1391f96d37fcSHong Zhang   merge->bj        = bj;
1392f96d37fcSHong Zhang   merge->coi       = coi;
1393f96d37fcSHong Zhang   merge->coj       = coj;
1394f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1395f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1396f96d37fcSHong Zhang   merge->owners_co = owners_co;
1397f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1398f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1399f96d37fcSHong Zhang 
14006da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1401c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1402*c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1403f96d37fcSHong Zhang 
1404f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1405f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1406f96d37fcSHong Zhang   c->ptap     = ptap;
14070298fd71SBarry Smith   ptap->api   = NULL;
14080298fd71SBarry Smith   ptap->apj   = NULL;
1409f96d37fcSHong Zhang   ptap->merge = merge;
1410d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1411d6ab1dc0SHong Zhang 
1412d6ab1dc0SHong Zhang   *C = Cmpi;
1413d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1414d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
141557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
141657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1417d6ab1dc0SHong Zhang   } else {
1418d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1419d6ab1dc0SHong Zhang   }
1420d6ab1dc0SHong Zhang #endif
1421d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1422d6ab1dc0SHong Zhang }
1423d6ab1dc0SHong Zhang 
1424d6ab1dc0SHong Zhang #undef __FUNCT__
14256da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14266da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1427d6ab1dc0SHong Zhang {
1428d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1429d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1430d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1431d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1432d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1433e745005bSHong Zhang   PetscInt            *adj;
1434e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1435e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1436d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1437ce94432eSBarry Smith   MPI_Comm            comm;
1438d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1439d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1440d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1441d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1442d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1443d6ab1dc0SHong Zhang   MPI_Status          *status;
1444d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1445d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1446a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1447d6ab1dc0SHong Zhang   Mat                 A_loc;
1448d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1449d6ab1dc0SHong Zhang 
1450d6ab1dc0SHong Zhang   PetscFunctionBegin;
1451ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1452d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1453d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1454d6ab1dc0SHong Zhang 
1455d6ab1dc0SHong Zhang   ptap  = c->ptap;
1456d6ab1dc0SHong Zhang   merge = ptap->merge;
1457d6ab1dc0SHong Zhang 
1458e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1459e745005bSHong Zhang   /*------------------------------------------*/
1460d6ab1dc0SHong Zhang   /* get data from symbolic products */
1461d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1462854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1463d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1464d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14651795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1466d6ab1dc0SHong Zhang 
1467d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1468d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1469d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1470d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1471d6ab1dc0SHong Zhang   ai    = a_loc->i;
1472d6ab1dc0SHong Zhang   aj    = a_loc->j;
1473d6ab1dc0SHong Zhang 
1474d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1475d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1476d6ab1dc0SHong Zhang     adj = aj + ai[i];
1477d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1478d6ab1dc0SHong Zhang 
1479d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1480e745005bSHong Zhang     /*-------------------------------------------------------------*/
1481d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1482d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1483d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1484d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1485d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1486d6ab1dc0SHong Zhang       row = poJ[j];
1487d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1488d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1489e745005bSHong Zhang       /* perform sparse axpy */
1490e745005bSHong Zhang       nexta  = 0;
1491d6ab1dc0SHong Zhang       valtmp = pA[j];
1492e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1493e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1494e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1495e745005bSHong Zhang           nexta++;
1496d6ab1dc0SHong Zhang         }
1497e745005bSHong Zhang       }
1498e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1499d6ab1dc0SHong Zhang     }
1500d6ab1dc0SHong Zhang 
1501d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1502d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1503d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1504d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1505d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1506d6ab1dc0SHong Zhang       row = pdJ[j];
1507d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1508d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1509e745005bSHong Zhang       /* perform sparse axpy */
1510e745005bSHong Zhang       nexta  = 0;
1511d6ab1dc0SHong Zhang       valtmp = pA[j];
1512e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1513e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1514e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1515e745005bSHong Zhang           nexta++;
1516d6ab1dc0SHong Zhang         }
1517e745005bSHong Zhang       }
1518e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1519d6ab1dc0SHong Zhang     }
1520d6ab1dc0SHong Zhang   }
1521d6ab1dc0SHong Zhang 
1522d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1523d6ab1dc0SHong Zhang   /*------------------------------------*/
1524d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1525d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1526d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1527d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1528d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang 
1530dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1531d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1532d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1533d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1534d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1535d6ab1dc0SHong Zhang     k++;
1536d6ab1dc0SHong Zhang   }
1537d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1538d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1539d6ab1dc0SHong Zhang 
1540d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1541d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1542d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1543d6ab1dc0SHong Zhang 
1544d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1545d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1546dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1547d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1548e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1549d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1550d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1551d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1552d6ab1dc0SHong Zhang   }
1553d6ab1dc0SHong Zhang 
1554d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1555d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1556d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1557d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1558d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1559d6ab1dc0SHong Zhang     /* add received vals into ba */
1560d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1561d6ab1dc0SHong Zhang       /* i-th row */
1562d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1563d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1564d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1565d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1566d6ab1dc0SHong Zhang         nextcj = 0;
1567d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1568d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1569d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1570d6ab1dc0SHong Zhang           }
1571d6ab1dc0SHong Zhang         }
1572d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1573d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1574d6ab1dc0SHong Zhang       }
1575d6ab1dc0SHong Zhang     }
1576d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1577d6ab1dc0SHong Zhang   }
1578d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1579d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1580d6ab1dc0SHong Zhang 
1581d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1582d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1583d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1584d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1585d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1586d6ab1dc0SHong Zhang }
1587d6ab1dc0SHong Zhang 
1588*c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
15892bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
15902bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1591d6ab1dc0SHong Zhang #undef __FUNCT__
15926da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
15936da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1594d6ab1dc0SHong Zhang {
1595d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1596d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1597d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15980298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1599d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1600d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1601d6ab1dc0SHong Zhang   PetscInt            nnz;
1602d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1603d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1604ce94432eSBarry Smith   MPI_Comm            comm;
1605d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1606d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1607d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1608d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1609d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1610d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1611d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1612d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1613d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1614d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1615c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1616d6ab1dc0SHong Zhang   PetscScalar         *vals;
1617d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1618d6ab1dc0SHong Zhang 
1619d6ab1dc0SHong Zhang   PetscFunctionBegin;
1620ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1621d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1622d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1623d6ab1dc0SHong 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);
1624d6ab1dc0SHong Zhang   }
1625d6ab1dc0SHong Zhang 
1626d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1627d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1628d6ab1dc0SHong Zhang 
1629d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1630b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1631d6ab1dc0SHong Zhang 
1632d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1633d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16342205254eSKarl Rupp 
1635d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1636d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1637d6ab1dc0SHong Zhang   ai          = a_loc->i;
1638d6ab1dc0SHong Zhang   aj          = a_loc->j;
1639d6ab1dc0SHong Zhang 
1640d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1641d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1642d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1643d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1644d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1645d6ab1dc0SHong Zhang 
1646d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1647d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1648d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1649d6ab1dc0SHong Zhang 
1650d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1651d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1652d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1653854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1654d6ab1dc0SHong Zhang   coi[0] = 0;
1655d6ab1dc0SHong Zhang 
1656d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1657d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1658a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1659d6ab1dc0SHong Zhang   current_space = free_space;
166019f0e57aSHong Zhang 
1661d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1662d6ab1dc0SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1663c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1664d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1665d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1666d6ab1dc0SHong Zhang 
1667d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1668d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1669d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1670d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1671d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1672d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1673d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1674d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1675d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1676d6ab1dc0SHong Zhang     }
1677d6ab1dc0SHong Zhang     nnz = lnk[0];
1678d6ab1dc0SHong Zhang 
1679d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1680d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1681d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1682d6ab1dc0SHong Zhang       nspacedouble++;
1683d6ab1dc0SHong Zhang     }
1684d6ab1dc0SHong Zhang 
1685d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1686d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16872205254eSKarl Rupp 
1688d6ab1dc0SHong Zhang     current_space->array           += nnz;
1689d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1690d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16912205254eSKarl Rupp 
1692d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1693d6ab1dc0SHong Zhang   }
1694d6ab1dc0SHong Zhang 
1695854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1696d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16972205254eSKarl Rupp 
1698118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1699d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1700d6ab1dc0SHong Zhang 
1701d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1702d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1703d6ab1dc0SHong Zhang   /* determine row ownership */
1704b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1705d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17062205254eSKarl Rupp 
1707d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1708d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17092205254eSKarl Rupp 
1710d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1711d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1712d6ab1dc0SHong Zhang 
1713d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17141795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1715785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17162205254eSKarl Rupp 
1717d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1718d6ab1dc0SHong Zhang   merge->nsend = 0;
1719d6ab1dc0SHong Zhang 
1720854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1721d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1722d6ab1dc0SHong Zhang 
1723d6ab1dc0SHong Zhang   proc = 0;
1724d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1725d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1726d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1727d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1728d6ab1dc0SHong Zhang   }
1729d6ab1dc0SHong Zhang 
1730d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1731d6ab1dc0SHong Zhang   owners_co[0] = 0;
1732d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1733d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1734d6ab1dc0SHong Zhang     if (len_si[proc]) {
1735d6ab1dc0SHong Zhang       merge->nsend++;
1736d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1737d6ab1dc0SHong Zhang       len         += len_si[proc];
1738d6ab1dc0SHong Zhang     }
1739d6ab1dc0SHong Zhang   }
1740d6ab1dc0SHong Zhang 
1741d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17420298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1743d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1744d6ab1dc0SHong Zhang 
1745d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1746d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1747d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1748854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1749d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1750d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1751d6ab1dc0SHong Zhang     i    = owners_co[proc];
1752d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1753d6ab1dc0SHong Zhang     k++;
1754d6ab1dc0SHong Zhang   }
1755d6ab1dc0SHong Zhang 
1756d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1757785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1758d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1759c280ed6aSJed Brown     PetscMPIInt icompleted;
1760c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1761d6ab1dc0SHong Zhang   }
1762d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1763d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1764d6ab1dc0SHong Zhang 
1765d6ab1dc0SHong Zhang   /* send and recv coi */
1766d6ab1dc0SHong Zhang   /*-------------------*/
1767d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1768d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1769854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1770d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1771d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1772d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1773d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1774d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1775d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1776d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1777d6ab1dc0SHong Zhang     */
1778d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1779d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1780d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1781d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1782d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1783d6ab1dc0SHong Zhang     nrows       = 0;
1784d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1785d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1786d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1787d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1788d6ab1dc0SHong Zhang       nrows++;
1789d6ab1dc0SHong Zhang     }
1790d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1791d6ab1dc0SHong Zhang     k++;
1792d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1793d6ab1dc0SHong Zhang   }
1794d6ab1dc0SHong Zhang   i = merge->nrecv;
1795d6ab1dc0SHong Zhang   while (i--) {
1796c280ed6aSJed Brown     PetscMPIInt icompleted;
1797c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1798d6ab1dc0SHong Zhang   }
1799d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1800d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1801d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1802d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1803d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1804d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1805d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1806d6ab1dc0SHong Zhang 
1807d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1808d6ab1dc0SHong Zhang   /*------------------------------------------*/
1809d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1810854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1811d6ab1dc0SHong Zhang   bi[0] = 0;
1812d6ab1dc0SHong Zhang 
1813d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1814d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1815a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1816d6ab1dc0SHong Zhang   current_space = free_space;
1817d6ab1dc0SHong Zhang 
1818dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1820d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1821d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1822d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18232205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1824d6ab1dc0SHong Zhang   }
1825d6ab1dc0SHong Zhang 
1826d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1827d6ab1dc0SHong Zhang   rmax = 0;
1828d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1829d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1830d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1831d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1832d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1833d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1834d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1835d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1836d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1837d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1838d6ab1dc0SHong Zhang     }
1839d6ab1dc0SHong Zhang 
1840d6ab1dc0SHong Zhang     /* add received col data into lnk */
1841d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1842d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1843d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1844d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1845d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1846d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1847d6ab1dc0SHong Zhang       }
1848d6ab1dc0SHong Zhang     }
1849d6ab1dc0SHong Zhang     nnz = lnk[0];
1850d6ab1dc0SHong Zhang 
1851d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1852d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1853d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1854d6ab1dc0SHong Zhang       nspacedouble++;
1855d6ab1dc0SHong Zhang     }
1856d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1857d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1858d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18592205254eSKarl Rupp 
1860d6ab1dc0SHong Zhang     current_space->array           += nnz;
1861d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1862d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18632205254eSKarl Rupp 
1864d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1865d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1866d6ab1dc0SHong Zhang   }
1867d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1868d6ab1dc0SHong Zhang 
1869854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1870d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1871118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1872d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1873d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1874d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1875d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1876d6ab1dc0SHong Zhang 
1877d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1878d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
18791795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1880d6ab1dc0SHong Zhang 
1881d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1882d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
188333d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1884d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1885d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1886d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1887d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1888d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1889d6ab1dc0SHong Zhang     row  = i + rstart;
1890d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1891d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1892d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1893d6ab1dc0SHong Zhang   }
1894d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1895d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1896d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1897d6ab1dc0SHong Zhang 
1898d6ab1dc0SHong Zhang   merge->bi        = bi;
1899d6ab1dc0SHong Zhang   merge->bj        = bj;
1900d6ab1dc0SHong Zhang   merge->coi       = coi;
1901d6ab1dc0SHong Zhang   merge->coj       = coj;
1902d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1903d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1904d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1905d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1906d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1907d6ab1dc0SHong Zhang 
19086da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1909d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1910*c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1911d6ab1dc0SHong Zhang 
1912d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1913d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19142205254eSKarl Rupp 
1915d6ab1dc0SHong Zhang   c->ptap     = ptap;
19160298fd71SBarry Smith   ptap->api   = NULL;
19170298fd71SBarry Smith   ptap->apj   = NULL;
1918d6ab1dc0SHong Zhang   ptap->merge = merge;
1919d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
19200298fd71SBarry Smith   ptap->apa   = NULL;
1921f96d37fcSHong Zhang 
1922f96d37fcSHong Zhang   *C = Cmpi;
1923f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1924f96d37fcSHong Zhang   if (bi[pn] != 0) {
192557622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
192657622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1927f96d37fcSHong Zhang   } else {
1928f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1929f96d37fcSHong Zhang   }
1930f96d37fcSHong Zhang #endif
1931f96d37fcSHong Zhang   PetscFunctionReturn(0);
1932f96d37fcSHong Zhang }
1933