xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 968382fd6276caa3ed3a1e5a3f88bb3b5c8dfc7f)
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;
121*968382fdSHong Zhang   } else {
122*968382fdSHong 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;
2469467f45fSHong Zhang   }
247598bc09dSHong Zhang 
248598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
249598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
250785e854fSJed Brown   ierr      = PetscMalloc1((am+2),&api);CHKERRQ(ierr);
251a1a4e84aSHong Zhang   ptap->api = api;
252598bc09dSHong Zhang   api[0]    = 0;
253598bc09dSHong Zhang 
254598bc09dSHong Zhang   /* create and initialize a linked list */
255f84c536bSHong Zhang   armax    = ad->rmax+ao->rmax;
2569467f45fSHong Zhang   if (size >1) {
257f84c536bSHong Zhang     prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
2589467f45fSHong Zhang   } else {
2599467f45fSHong Zhang     prmax = p_loc->rmax;
2609467f45fSHong Zhang   }
261f84c536bSHong Zhang   nlnk_max = armax*prmax;
262f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2630d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
264598bc09dSHong Zhang 
265bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
266bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
2672205254eSKarl Rupp 
268598bc09dSHong Zhang   current_space = free_space;
269598bc09dSHong Zhang 
270598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
271598bc09dSHong Zhang   for (i=0; i<am; i++) {
272598bc09dSHong Zhang     /* diagonal portion of A */
273598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
274598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
275598bc09dSHong Zhang       row  = *adj++;
276598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
277598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
278598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2791fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
280598bc09dSHong Zhang     }
281598bc09dSHong Zhang     /* off-diagonal portion of A */
282598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
283598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
284598bc09dSHong Zhang       row  = *aoj++;
285598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
286598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2871fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
288598bc09dSHong Zhang     }
289598bc09dSHong Zhang 
290d14fa243SHong Zhang     apnz     = lnk[0];
291598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
292598bc09dSHong Zhang 
293598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
294598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
295598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
296598bc09dSHong Zhang       nspacedouble++;
297598bc09dSHong Zhang     }
298598bc09dSHong Zhang 
299598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
3001fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
301598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
3022205254eSKarl Rupp 
303598bc09dSHong Zhang     current_space->array           += apnz;
304598bc09dSHong Zhang     current_space->local_used      += apnz;
305598bc09dSHong Zhang     current_space->local_remaining -= apnz;
306598bc09dSHong Zhang   }
307598bc09dSHong Zhang 
308598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
309598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
310785e854fSJed Brown   ierr = PetscMalloc1((api[am]+1),&ptap->apj);CHKERRQ(ierr);
311a1a4e84aSHong Zhang   apj  = ptap->apj;
312a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
313598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
314598bc09dSHong Zhang 
315f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
3161795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
3172205254eSKarl Rupp 
318f84c536bSHong Zhang   ptap->apa = apa;
319f84c536bSHong Zhang 
320bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
321598bc09dSHong Zhang   /*----------------------------------------------------*/
322bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
323bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
32433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
325a2f3521dSMark F. Adams 
326bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
327bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
328598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
329598bc09dSHong Zhang   for (i=0; i<am; i++) {
330598bc09dSHong Zhang     row  = i + rstart;
331598bc09dSHong Zhang     apnz = api[i+1] - api[i];
332bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
333598bc09dSHong Zhang     apj += apnz;
334598bc09dSHong Zhang   }
335bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
336bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
337598bc09dSHong Zhang 
338bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
339bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
340f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
341bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
342bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
343598bc09dSHong Zhang 
344bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
345bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
346598bc09dSHong Zhang   c->ptap = ptap;
347598bc09dSHong Zhang 
348bfcd1a73SHong Zhang   *C = Cmpi;
349bfcd1a73SHong Zhang 
350bfcd1a73SHong Zhang   /* set MatInfo */
351118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
352bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
353bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
354bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
355bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
356bfcd1a73SHong Zhang 
357bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
358bfcd1a73SHong Zhang   if (api[am]) {
35957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
36057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
361bfcd1a73SHong Zhang   } else {
362bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
363bfcd1a73SHong Zhang   }
364bfcd1a73SHong Zhang #endif
365598bc09dSHong Zhang   PetscFunctionReturn(0);
366598bc09dSHong Zhang }
367598bc09dSHong Zhang 
3689123193aSHong Zhang #undef __FUNCT__
3699123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3709123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3719123193aSHong Zhang {
3729123193aSHong Zhang   PetscErrorCode ierr;
3739123193aSHong Zhang 
3749123193aSHong Zhang   PetscFunctionBegin;
3759123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3763ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3779123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3783ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3799123193aSHong Zhang   }
3803ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3819123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3823ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3839123193aSHong Zhang   PetscFunctionReturn(0);
3849123193aSHong Zhang }
3859123193aSHong Zhang 
3864324174eSBarry Smith typedef struct {
3874324174eSBarry Smith   Mat         workB;
3884324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3894324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3904324174eSBarry Smith } MPIAIJ_MPIDense;
3914324174eSBarry Smith 
3927af9e4fdSHong Zhang #undef __FUNCT__
39396b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
39496b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3954324174eSBarry Smith {
3964324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3974324174eSBarry Smith   PetscErrorCode  ierr;
3984324174eSBarry Smith 
3994324174eSBarry Smith   PetscFunctionBegin;
4006bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
4014324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
4024324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
4034324174eSBarry Smith   PetscFunctionReturn(0);
4044324174eSBarry Smith }
4054324174eSBarry Smith 
4069123193aSHong Zhang #undef __FUNCT__
4079123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4089123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4099123193aSHong Zhang {
4109123193aSHong Zhang   PetscErrorCode         ierr;
411f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
412d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
413bf0cc555SLisandro Dalcin   PetscContainer         container;
4144324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4154324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4164324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4174324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
418d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4199123193aSHong Zhang 
4209123193aSHong Zhang   PetscFunctionBegin;
421ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
422d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
42333d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
424cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4250298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
426cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
427cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4282205254eSKarl Rupp 
429d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
430f32d5d43SBarry Smith 
431b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
432f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4330298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4344324174eSBarry Smith   /* Create work arrays needed */
435dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
436dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
437dcca6d9dSJed Brown                       from->n,&contents->rwaits,
438dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4394324174eSBarry Smith 
440ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
441bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
44296b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
443bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
444bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4459123193aSHong Zhang   PetscFunctionReturn(0);
4469123193aSHong Zhang }
4479123193aSHong Zhang 
4487af9e4fdSHong Zhang #undef __FUNCT__
4497af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
450f32d5d43SBarry Smith /*
4512636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4522636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
453f32d5d43SBarry Smith */
4544324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
455f32d5d43SBarry Smith {
456f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
457f32d5d43SBarry Smith   PetscErrorCode         ierr;
458f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
459f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
460f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
461f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
462f32d5d43SBarry Smith   PetscInt               i,j,k;
463aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
464aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
465f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
466ce94432eSBarry Smith   MPI_Comm               comm;
467d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
468f32d5d43SBarry Smith   MPI_Status             status;
4694324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
470bf0cc555SLisandro Dalcin   PetscContainer         container;
4714324174eSBarry Smith   Mat                    workB;
472f32d5d43SBarry Smith 
473f32d5d43SBarry Smith   PetscFunctionBegin;
474ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
475bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
47629825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
477bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4784324174eSBarry Smith 
4794324174eSBarry Smith   workB = *outworkB = contents->workB;
480cf1a0672SHong 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);
481f32d5d43SBarry Smith   sindices = to->indices;
482f32d5d43SBarry Smith   sstarts  = to->starts;
483f32d5d43SBarry Smith   sprocs   = to->procs;
4844324174eSBarry Smith   swaits   = contents->swaits;
4854324174eSBarry Smith   svalues  = contents->svalues;
486f32d5d43SBarry Smith 
487f32d5d43SBarry Smith   rindices = from->indices;
488f32d5d43SBarry Smith   rstarts  = from->starts;
489f32d5d43SBarry Smith   rprocs   = from->procs;
4904324174eSBarry Smith   rwaits   = contents->rwaits;
4914324174eSBarry Smith   rvalues  = contents->rvalues;
492f32d5d43SBarry Smith 
4938c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
4948c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
495f32d5d43SBarry Smith 
496f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
497f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
498f32d5d43SBarry Smith   }
499f32d5d43SBarry Smith 
500f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
501f32d5d43SBarry Smith     /* pack a message at a time */
502f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
503f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
504f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5052636ff26SBarry Smith       }
5062636ff26SBarry Smith     }
507f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
508f32d5d43SBarry Smith   }
5092636ff26SBarry Smith 
510f32d5d43SBarry Smith   nrecvs = from->n;
511f32d5d43SBarry Smith   while (nrecvs) {
512f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
513f32d5d43SBarry Smith     nrecvs--;
514f32d5d43SBarry Smith     /* unpack a message at a time */
515f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
516f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
517f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5182636ff26SBarry Smith       }
5192636ff26SBarry Smith     }
520f32d5d43SBarry Smith   }
521cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
522f32d5d43SBarry Smith 
5238c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5248c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
525f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
526f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
527f32d5d43SBarry Smith   PetscFunctionReturn(0);
528f32d5d43SBarry Smith }
529f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
530f32d5d43SBarry Smith 
5319123193aSHong Zhang #undef __FUNCT__
5329123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5339123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5349123193aSHong Zhang {
5359123193aSHong Zhang   PetscErrorCode ierr;
536f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
537f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
538f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
539f32d5d43SBarry Smith   Mat            workB;
5409123193aSHong Zhang 
5419123193aSHong Zhang   PetscFunctionBegin;
542f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
543f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
544f32d5d43SBarry Smith 
545f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5464324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
547f32d5d43SBarry Smith 
548f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
549f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5509123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5519123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5529123193aSHong Zhang   PetscFunctionReturn(0);
5539123193aSHong Zhang }
554cf1a0672SHong Zhang 
5551634b032SHong Zhang #undef __FUNCT__
556f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
557f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5581634b032SHong Zhang {
5591634b032SHong Zhang   PetscErrorCode ierr;
560cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
561cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
562cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
563cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
564cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
565cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
566cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
56729825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
56829825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
569cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
57029825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
571cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
57229825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
57329825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
574a7c7454dSHong Zhang   MPI_Comm       comm;
575a7c7454dSHong Zhang   PetscMPIInt    size;
5761634b032SHong Zhang 
5771634b032SHong Zhang   PetscFunctionBegin;
578a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
579a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
580a7c7454dSHong Zhang 
581cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
582cf1a0672SHong Zhang   /*-----------------------------------------------------*/
583cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
584b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
585cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
586cf1a0672SHong Zhang 
587cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
588cf1a0672SHong Zhang   /*----------------------------------------------------------*/
589cf1a0672SHong Zhang   /* get data from symbolic products */
590cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
591cf1a0672SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
592a7c7454dSHong Zhang   if (size >1) {
593a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
594cf1a0672SHong Zhang     pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
595a7c7454dSHong Zhang   }
596cf1a0672SHong Zhang 
59729825010SHong Zhang   api = ptap->api;
59829825010SHong Zhang   apj = ptap->apj;
599cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
60029825010SHong Zhang     apJ = apj + api[i];
60129825010SHong Zhang 
602cf1a0672SHong Zhang     /* diagonal portion of A */
603cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
604cf1a0672SHong Zhang     adj = ad->j + adi[i];
605cf1a0672SHong Zhang     ada = ad->a + adi[i];
606cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
607cf1a0672SHong Zhang       row = adj[j];
608cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
609cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
610cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
61129825010SHong Zhang       /* perform sparse axpy */
612cf1a0672SHong Zhang       valtmp = ada[j];
61329825010SHong Zhang       nextp  = 0;
61429825010SHong Zhang       for (k=0; nextp<pnz; k++) {
61529825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
61629825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
61729825010SHong Zhang         }
6181634b032SHong Zhang       }
619cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
620cf1a0672SHong Zhang     }
6211634b032SHong Zhang 
622cf1a0672SHong Zhang     /* off-diagonal portion of A */
623cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
624cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
625cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
626cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
627cf1a0672SHong Zhang       row = aoj[j];
628cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
629cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
630cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
63129825010SHong Zhang       /* perform sparse axpy */
632cf1a0672SHong Zhang       valtmp = aoa[j];
63329825010SHong Zhang       nextp  = 0;
63429825010SHong Zhang       for (k=0; nextp<pnz; k++) {
63529825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
63629825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
63729825010SHong Zhang         }
638cf1a0672SHong Zhang       }
639cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
640cf1a0672SHong Zhang     }
6411634b032SHong Zhang 
642cf1a0672SHong Zhang     /* set values in C */
643cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
644cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6451634b032SHong Zhang 
646cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
647cf1a0672SHong Zhang     ca = coa + co->i[i];
648cf1a0672SHong Zhang     k  = 0;
649cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
650cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
65129825010SHong Zhang       ca[k0]        = apa_sparse[k];
65229825010SHong Zhang       apa_sparse[k] = 0.0;
653cf1a0672SHong Zhang       k++;
654cf1a0672SHong Zhang     }
6551634b032SHong Zhang 
656cf1a0672SHong Zhang     /* diagonal part of C */
657cf1a0672SHong Zhang     ca = cda + cd->i[i];
658cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
65929825010SHong Zhang       ca[k1]        = apa_sparse[k];
66029825010SHong Zhang       apa_sparse[k] = 0.0;
661cf1a0672SHong Zhang       k++;
662cf1a0672SHong Zhang     }
663cf1a0672SHong Zhang 
664cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
665cf1a0672SHong Zhang     ca = coa + co->i[i];
666cf1a0672SHong Zhang     for (; k0<conz; k0++) {
66729825010SHong Zhang       ca[k0]        = apa_sparse[k];
66829825010SHong Zhang       apa_sparse[k] = 0.0;
669cf1a0672SHong Zhang       k++;
670cf1a0672SHong Zhang     }
671cf1a0672SHong Zhang   }
672cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
673cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
674cf1a0672SHong Zhang   PetscFunctionReturn(0);
675cf1a0672SHong Zhang }
676cf1a0672SHong Zhang 
6770fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
678cf1a0672SHong Zhang #undef __FUNCT__
679b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
680b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
681cf1a0672SHong Zhang {
682cf1a0672SHong Zhang   PetscErrorCode     ierr;
683ce94432eSBarry Smith   MPI_Comm           comm;
684a7c7454dSHong Zhang   PetscMPIInt        size;
685cf1a0672SHong Zhang   Mat                Cmpi;
686cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6870298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
688cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
689cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
690cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
691cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
692f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
693cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
694cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
695cf1a0672SHong Zhang   PetscReal          afill;
696cf1a0672SHong Zhang   PetscScalar        *apa;
697cf1a0672SHong Zhang 
698cf1a0672SHong Zhang   PetscFunctionBegin;
699ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
700a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
701a7c7454dSHong Zhang 
702cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
703b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
704cf1a0672SHong Zhang 
705cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
706b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
70719f0e57aSHong Zhang 
708cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
709cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
710cf1a0672SHong Zhang 
711cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
712cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
713a7c7454dSHong Zhang   if (size > 1) {
714a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
715a7c7454dSHong Zhang     pi_oth = p_oth->i;
716a7c7454dSHong Zhang     pj_oth = p_oth->j;
717*968382fdSHong Zhang   } else {
718*968382fdSHong Zhang     pi_oth = NULL;
719*968382fdSHong Zhang     pj_oth = NULL;
720a7c7454dSHong Zhang   }
721cf1a0672SHong Zhang 
722cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
723cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
724785e854fSJed Brown   ierr      = PetscMalloc1((am+2),&api);CHKERRQ(ierr);
725cf1a0672SHong Zhang   ptap->api = api;
726cf1a0672SHong Zhang   api[0]    = 0;
727cf1a0672SHong Zhang 
728cf1a0672SHong Zhang   /* create and initialize a linked list */
729cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
730a7c7454dSHong Zhang   if (size >1) {
731cf1a0672SHong Zhang     prmax = PetscMax(p_loc->rmax,p_oth->rmax);
732a7c7454dSHong Zhang   } else {
733a7c7454dSHong Zhang     prmax = p_loc->rmax;
734a7c7454dSHong Zhang   }
735cf1a0672SHong Zhang   nlnk_max = armax*prmax;
736cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
737f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
738cf1a0672SHong Zhang 
739cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
740cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
7412205254eSKarl Rupp 
742cf1a0672SHong Zhang   current_space = free_space;
743cf1a0672SHong Zhang 
744cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
745cf1a0672SHong Zhang   for (i=0; i<am; i++) {
746cf1a0672SHong Zhang     /* diagonal portion of A */
747cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
748cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
749cf1a0672SHong Zhang       row  = *adj++;
750cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
751cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
752cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
753f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
754cf1a0672SHong Zhang     }
755cf1a0672SHong Zhang     /* off-diagonal portion of A */
756cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
757cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
758cf1a0672SHong Zhang       row  = *aoj++;
759cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
760cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
761f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
762cf1a0672SHong Zhang     }
763cf1a0672SHong Zhang 
764f84c536bSHong Zhang     apnz     = *lnk;
765cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
766cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
767cf1a0672SHong Zhang 
768cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
769cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
770cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
771cf1a0672SHong Zhang       nspacedouble++;
772cf1a0672SHong Zhang     }
773cf1a0672SHong Zhang 
774cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
775f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
776cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7772205254eSKarl Rupp 
778cf1a0672SHong Zhang     current_space->array           += apnz;
779cf1a0672SHong Zhang     current_space->local_used      += apnz;
780cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
781cf1a0672SHong Zhang   }
782cf1a0672SHong Zhang 
783cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
784cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
785785e854fSJed Brown   ierr = PetscMalloc1((api[am]+1),&ptap->apj);CHKERRQ(ierr);
786cf1a0672SHong Zhang   apj  = ptap->apj;
787cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
788f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
789cf1a0672SHong Zhang 
790cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
791cf1a0672SHong Zhang   /*----------------------------------------------------*/
792cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
793cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
79433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
795cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
796cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
797cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
798cf1a0672SHong Zhang 
79929825010SHong Zhang   /* malloc apa for assembly Cmpi */
8001795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8012205254eSKarl Rupp 
80229825010SHong Zhang   ptap->apa = apa;
803cf1a0672SHong Zhang   for (i=0; i<am; i++) {
804cf1a0672SHong Zhang     row  = i + rstart;
805cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
806cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
807cf1a0672SHong Zhang     apj += apnz;
808cf1a0672SHong Zhang   }
809cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
810cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
811cf1a0672SHong Zhang 
812cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
813cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
814f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
815cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
816cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
817cf1a0672SHong Zhang 
818cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
819cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
820cf1a0672SHong Zhang   c->ptap = ptap;
821cf1a0672SHong Zhang 
822cf1a0672SHong Zhang   *C = Cmpi;
823cf1a0672SHong Zhang 
824cf1a0672SHong Zhang   /* set MatInfo */
825118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
826cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
827cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
828cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
829cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
830cf1a0672SHong Zhang 
831cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
832cf1a0672SHong Zhang   if (api[am]) {
83357622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
83457622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
835cf1a0672SHong Zhang   } else {
836cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
837cf1a0672SHong Zhang   }
838cf1a0672SHong Zhang #endif
8391634b032SHong Zhang   PetscFunctionReturn(0);
8401634b032SHong Zhang }
841f96d37fcSHong Zhang 
842f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
843f96d37fcSHong Zhang #undef __FUNCT__
844f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
845c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
846f96d37fcSHong Zhang {
847f96d37fcSHong Zhang   PetscErrorCode ierr;
848c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
849c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
850f96d37fcSHong Zhang 
851f96d37fcSHong Zhang   PetscFunctionBegin;
852c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
853c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
854c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
855c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
856c216dbf3SHong Zhang 
857c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
858c216dbf3SHong Zhang     switch (alg) {
859c216dbf3SHong Zhang     case 1:
8602bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
861c216dbf3SHong Zhang       break;
862c216dbf3SHong Zhang     case 2:
863275476c6SMatthew G. Knepley     {
864ab1b013aSHong Zhang       Mat         Pt;
865ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
866ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
867ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
868ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
869ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
870ab1b013aSHong Zhang       ptap     = c->ptap;
871ab1b013aSHong Zhang       ptap->Pt = Pt;
872c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
873c216dbf3SHong Zhang       PetscFunctionReturn(0);
874275476c6SMatthew G. Knepley     }
875c216dbf3SHong Zhang       break;
876c216dbf3SHong Zhang     default:
8776da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
878c216dbf3SHong Zhang       break;
879c216dbf3SHong Zhang     }
880c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
881c216dbf3SHong Zhang   }
882c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
883c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
884c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
885ab1b013aSHong Zhang   PetscFunctionReturn(0);
886ab1b013aSHong Zhang }
887ab1b013aSHong Zhang 
888c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
889c216dbf3SHong Zhang #undef __FUNCT__
890c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
891c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
892c216dbf3SHong Zhang {
893c216dbf3SHong Zhang   PetscErrorCode ierr;
8942bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
8952bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
8962bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
897c216dbf3SHong Zhang 
898c216dbf3SHong Zhang   PetscFunctionBegin;
899c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
900c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
901f96d37fcSHong Zhang   PetscFunctionReturn(0);
902f96d37fcSHong Zhang }
903f96d37fcSHong Zhang 
9046da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
905f96d37fcSHong Zhang #undef __FUNCT__
9066da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9076da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
908f96d37fcSHong Zhang {
909c5af039cSHong Zhang   PetscErrorCode      ierr;
910c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
911497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
912c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
913c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
914d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
915497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
916e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
917c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
918ce94432eSBarry Smith   MPI_Comm            comm;
919c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
920c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
921c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
922c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
923c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
924c5af039cSHong Zhang   MPI_Status          *status;
925c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
926d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
927a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
928c5af039cSHong Zhang   Mat                 A_loc;
929c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
930f96d37fcSHong Zhang 
931f96d37fcSHong Zhang   PetscFunctionBegin;
932ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
933c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
934c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
935c5af039cSHong Zhang 
936c5af039cSHong Zhang   ptap  = c->ptap;
937c5af039cSHong Zhang   merge = ptap->merge;
938c5af039cSHong Zhang 
939c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
940c5af039cSHong Zhang   /*--------------------------------------------------------------*/
941c5af039cSHong Zhang   /* get data from symbolic products */
942c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
9431795a4d1SJed Brown   ierr = PetscCalloc1((coi[pon]+1),&coa);CHKERRQ(ierr);
944c5af039cSHong Zhang 
945c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
946c5af039cSHong Zhang   owners = merge->rowmap->range;
9471795a4d1SJed Brown   ierr   = PetscCalloc1((bi[cm]+1),&ba);CHKERRQ(ierr);
948c5af039cSHong Zhang 
949c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
950c5af039cSHong Zhang   A_loc = ptap->A_loc;
951c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
952c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
953d6ab1dc0SHong Zhang   ai    = a_loc->i;
954d6ab1dc0SHong Zhang   aj    = a_loc->j;
955c5af039cSHong Zhang 
9561795a4d1SJed Brown   ierr = PetscCalloc1((A->cmap->N),&aval);CHKERRQ(ierr); /* non-scalable!!! */
957c5af039cSHong Zhang 
958c5af039cSHong Zhang   for (i=0; i<am; i++) {
959e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
960d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
961d6ab1dc0SHong Zhang     adj = aj + ai[i];
962d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
963c5af039cSHong Zhang     for (j=0; j<anz; j++) {
964e745005bSHong Zhang       aval[adj[j]] = ada[j];
965c5af039cSHong Zhang     }
966c5af039cSHong Zhang 
967c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
968c5af039cSHong Zhang     /*--------------------------------------------------------------*/
969c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
970c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
971c5af039cSHong Zhang     poJ = po->j + po->i[i];
972c5af039cSHong Zhang     pA  = po->a + po->i[i];
973c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
974c5af039cSHong Zhang       row = poJ[j];
975c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
976c5af039cSHong Zhang       cj  = coj + coi[row];
977c5af039cSHong Zhang       ca  = coa + coi[row];
978c5af039cSHong Zhang       /* perform dense axpy */
979c5af039cSHong Zhang       valtmp = pA[j];
980c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
981e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
982c5af039cSHong Zhang       }
983c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
984c5af039cSHong Zhang     }
985c5af039cSHong Zhang 
986c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
987c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
988c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
989c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
990c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
991c5af039cSHong Zhang       row = pdJ[j];
992c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
993c5af039cSHong Zhang       cj  = bj + bi[row];
994c5af039cSHong Zhang       ca  = ba + bi[row];
995c5af039cSHong Zhang       /* perform dense axpy */
996c5af039cSHong Zhang       valtmp = pA[j];
997c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
998e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
999c5af039cSHong Zhang       }
1000c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1001c5af039cSHong Zhang     }
1002c5af039cSHong Zhang 
1003d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1004d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1005e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1006c5af039cSHong Zhang   }
1007c5af039cSHong Zhang 
1008c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1009c5af039cSHong Zhang   /*------------------------------------*/
1010c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1011c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1012c5af039cSHong Zhang   len_s  = merge->len_s;
1013c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1014c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1015c5af039cSHong Zhang 
1016dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1017c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1018c5af039cSHong Zhang     if (!len_s[proc]) continue;
1019c5af039cSHong Zhang     i    = merge->owners_co[proc];
1020c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1021c5af039cSHong Zhang     k++;
1022c5af039cSHong Zhang   }
1023c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1024c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1025c5af039cSHong Zhang 
1026c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1027c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1028c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1029c5af039cSHong Zhang 
1030c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1031c5af039cSHong Zhang   /*----------------------------------------------------*/
1032dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1033c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1034c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1035c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1036c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1037c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1038c5af039cSHong Zhang   }
1039c5af039cSHong Zhang 
1040c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1041c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1042c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1043c5af039cSHong Zhang     ba_i = ba + bi[i];
1044c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1045c5af039cSHong Zhang     /* add received vals into ba */
1046c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1047c5af039cSHong Zhang       /* i-th row */
1048c5af039cSHong Zhang       if (i == *nextrow[k]) {
1049c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1050c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1051c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1052c5af039cSHong Zhang         nextcj = 0;
1053c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1054c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1055c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1056c5af039cSHong Zhang           }
1057c5af039cSHong Zhang         }
1058c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1059c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1060c5af039cSHong Zhang       }
1061c5af039cSHong Zhang     }
1062c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1063c5af039cSHong Zhang   }
1064c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1065c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1066c5af039cSHong Zhang 
1067c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1068c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1069c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1070c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1071e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1072f96d37fcSHong Zhang   PetscFunctionReturn(0);
1073f96d37fcSHong Zhang }
1074f96d37fcSHong Zhang 
1075c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1076f96d37fcSHong Zhang #undef __FUNCT__
10772bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
10782bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1079f96d37fcSHong Zhang {
1080f96d37fcSHong Zhang   PetscErrorCode      ierr;
10814e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1082f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10830298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1084f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1085f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1086f96d37fcSHong Zhang   PetscInt            nnz;
10874e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1088497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1089f96d37fcSHong Zhang   PetscBT             lnkbt;
1090ce94432eSBarry Smith   MPI_Comm            comm;
1091f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1092f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1093f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1094f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1095f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1096f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1097f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1098f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1099d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1100f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1101c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1102f96d37fcSHong Zhang   PetscScalar         *vals;
11034e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1104f96d37fcSHong Zhang 
1105f96d37fcSHong Zhang   PetscFunctionBegin;
1106ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1107c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1108c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1109c5af039cSHong 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);
1110c5af039cSHong Zhang   }
1111c5af039cSHong Zhang 
1112f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1113f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1114f96d37fcSHong Zhang 
1115f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1116b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1117f96d37fcSHong Zhang 
1118f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1119f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11202205254eSKarl Rupp 
1121c5af039cSHong Zhang   ptap->A_loc = A_loc;
11222205254eSKarl Rupp 
11231c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1124d6ab1dc0SHong Zhang   ai    = a_loc->i;
1125d6ab1dc0SHong Zhang   aj    = a_loc->j;
1126f96d37fcSHong Zhang 
1127f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1128f96d37fcSHong Zhang   /*----------------------------------------------------*/
11294e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11304e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11314e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11324e938277SHong Zhang 
11334e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11344e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11354e938277SHong Zhang   poti = pot->i; potj = pot->j;
1136f96d37fcSHong Zhang 
1137f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11382205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1139785e854fSJed Brown   ierr   = PetscMalloc1((pon+1),&coi);CHKERRQ(ierr);
1140f96d37fcSHong Zhang   coi[0] = 0;
1141f96d37fcSHong Zhang 
1142f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1143d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1144a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1145f96d37fcSHong Zhang   current_space = free_space;
114619f0e57aSHong Zhang 
1147f96d37fcSHong Zhang   /* create and initialize a linked list */
11484e938277SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1149c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
11504e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
11514e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1152f96d37fcSHong Zhang 
1153f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1154f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1155f96d37fcSHong Zhang     ptJ = potj + poti[i];
1156f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1157f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1158d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1159d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1160f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1161d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1162f96d37fcSHong Zhang     }
11634e938277SHong Zhang     nnz = lnk[0];
1164f96d37fcSHong Zhang 
1165f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1166f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1167f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1168f96d37fcSHong Zhang       nspacedouble++;
1169f96d37fcSHong Zhang     }
1170f96d37fcSHong Zhang 
1171f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11724e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11732205254eSKarl Rupp 
1174f96d37fcSHong Zhang     current_space->array           += nnz;
1175f96d37fcSHong Zhang     current_space->local_used      += nnz;
1176f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11772205254eSKarl Rupp 
1178f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1179f96d37fcSHong Zhang   }
1180f96d37fcSHong Zhang 
1181785e854fSJed Brown   ierr = PetscMalloc1((coi[pon]+1),&coj);CHKERRQ(ierr);
1182f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11832205254eSKarl Rupp 
1184118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1185f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1186f96d37fcSHong Zhang 
1187f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1188f96d37fcSHong Zhang   /*----------------------------------------------*/
1189f96d37fcSHong Zhang   /* determine row ownership */
1190b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1191f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11922205254eSKarl Rupp 
1193f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1194f96d37fcSHong Zhang   merge->rowmap->bs = 1;
11952205254eSKarl Rupp 
1196f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1197f96d37fcSHong Zhang   owners = merge->rowmap->range;
1198f96d37fcSHong Zhang 
1199f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12001795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1201785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12022205254eSKarl Rupp 
1203f96d37fcSHong Zhang   len_s        = merge->len_s;
1204f96d37fcSHong Zhang   merge->nsend = 0;
1205f96d37fcSHong Zhang 
1206785e854fSJed Brown   ierr = PetscMalloc1((size+2),&owners_co);CHKERRQ(ierr);
1207f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1208f96d37fcSHong Zhang 
1209f96d37fcSHong Zhang   proc = 0;
1210f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1211f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1212f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1213f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1214f96d37fcSHong Zhang   }
1215f96d37fcSHong Zhang 
1216f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1217f96d37fcSHong Zhang   owners_co[0] = 0;
1218f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1219f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1220f96d37fcSHong Zhang     if (len_si[proc]) {
1221f96d37fcSHong Zhang       merge->nsend++;
1222f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1223f96d37fcSHong Zhang       len         += len_si[proc];
1224f96d37fcSHong Zhang     }
1225f96d37fcSHong Zhang   }
1226f96d37fcSHong Zhang 
1227f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12280298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1229f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1230f96d37fcSHong Zhang 
1231f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1232f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1233f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1234785e854fSJed Brown   ierr = PetscMalloc1((merge->nsend+1),&swaits);CHKERRQ(ierr);
1235f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1236f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1237f96d37fcSHong Zhang     i    = owners_co[proc];
1238f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1239f96d37fcSHong Zhang     k++;
1240f96d37fcSHong Zhang   }
1241f96d37fcSHong Zhang 
1242f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1243785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1244f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1245c280ed6aSJed Brown     PetscMPIInt icompleted;
1246c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1247f96d37fcSHong Zhang   }
1248f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1249f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1250f96d37fcSHong Zhang 
1251f96d37fcSHong Zhang   /* send and recv coi */
1252f96d37fcSHong Zhang   /*-------------------*/
1253f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1254f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1255785e854fSJed Brown   ierr   = PetscMalloc1((len+1),&buf_s);CHKERRQ(ierr);
1256f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1257f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1258f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1259f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1260f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1261f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1262f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1263f96d37fcSHong Zhang     */
1264f96d37fcSHong Zhang     /*-------------------------------------------*/
1265f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1266f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1267f96d37fcSHong Zhang     buf_si[0]   = nrows;
1268f96d37fcSHong Zhang     buf_si_i[0] = 0;
1269f96d37fcSHong Zhang     nrows       = 0;
1270f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1271f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1272f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1273f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1274f96d37fcSHong Zhang       nrows++;
1275f96d37fcSHong Zhang     }
1276f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1277f96d37fcSHong Zhang     k++;
1278f96d37fcSHong Zhang     buf_si += len_si[proc];
1279f96d37fcSHong Zhang   }
1280f96d37fcSHong Zhang   i = merge->nrecv;
1281f96d37fcSHong Zhang   while (i--) {
1282c280ed6aSJed Brown     PetscMPIInt icompleted;
1283c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1284f96d37fcSHong Zhang   }
1285f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1286f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1287f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1288f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1289f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1290f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1291f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1292f96d37fcSHong Zhang 
1293f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1294f96d37fcSHong Zhang   /*------------------------------------------*/
1295f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1296785e854fSJed Brown   ierr  = PetscMalloc1((pn+1),&bi);CHKERRQ(ierr);
1297f96d37fcSHong Zhang   bi[0] = 0;
1298f96d37fcSHong Zhang 
1299c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1300d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1301a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1302f96d37fcSHong Zhang   current_space = free_space;
1303f96d37fcSHong Zhang 
1304dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1305f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1306f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1307f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1308f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1309f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1310f96d37fcSHong Zhang   }
1311f96d37fcSHong Zhang 
13121c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1313f96d37fcSHong Zhang   rmax = 0;
1314f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1315f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1316f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1317f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1318f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1319f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1320d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1321d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1322f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1323d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1324f96d37fcSHong Zhang     }
13254e938277SHong Zhang 
1326f96d37fcSHong Zhang     /* add received col data into lnk */
1327f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1328f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1329f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1330f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13314e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1332f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1333f96d37fcSHong Zhang       }
1334f96d37fcSHong Zhang     }
13354e938277SHong Zhang     nnz = lnk[0];
1336f96d37fcSHong Zhang 
1337f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1338f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1339f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1340f96d37fcSHong Zhang       nspacedouble++;
1341f96d37fcSHong Zhang     }
1342f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13434e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1344f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13452205254eSKarl Rupp 
1346f96d37fcSHong Zhang     current_space->array           += nnz;
1347f96d37fcSHong Zhang     current_space->local_used      += nnz;
1348f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13492205254eSKarl Rupp 
1350f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1351f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1352f96d37fcSHong Zhang   }
1353f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1354f96d37fcSHong Zhang 
1355785e854fSJed Brown   ierr = PetscMalloc1((bi[pn]+1),&bj);CHKERRQ(ierr);
1356f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13572205254eSKarl Rupp 
1358118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1359f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1360d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13614e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13624e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1363f96d37fcSHong Zhang 
13641c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13651c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13661795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1367f96d37fcSHong Zhang 
1368f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1369f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
137033d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1371f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1372f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1373f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1374f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1375f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1376f96d37fcSHong Zhang     row  = i + rstart;
1377f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1378f96d37fcSHong Zhang     Jptr = bj + bi[i];
1379f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1380f96d37fcSHong Zhang   }
1381f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1382f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1383f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1384f96d37fcSHong Zhang 
1385f96d37fcSHong Zhang   merge->bi        = bi;
1386f96d37fcSHong Zhang   merge->bj        = bj;
1387f96d37fcSHong Zhang   merge->coi       = coi;
1388f96d37fcSHong Zhang   merge->coj       = coj;
1389f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1390f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1391f96d37fcSHong Zhang   merge->owners_co = owners_co;
1392f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1393f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1394f96d37fcSHong Zhang 
13956da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1396c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1397f96d37fcSHong Zhang 
1398f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1399f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1400f96d37fcSHong Zhang   c->ptap     = ptap;
14010298fd71SBarry Smith   ptap->api   = NULL;
14020298fd71SBarry Smith   ptap->apj   = NULL;
1403f96d37fcSHong Zhang   ptap->merge = merge;
1404d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1405d6ab1dc0SHong Zhang 
1406d6ab1dc0SHong Zhang   *C = Cmpi;
1407d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1408d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
140957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
141057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1411d6ab1dc0SHong Zhang   } else {
1412d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1413d6ab1dc0SHong Zhang   }
1414d6ab1dc0SHong Zhang #endif
1415d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1416d6ab1dc0SHong Zhang }
1417d6ab1dc0SHong Zhang 
1418d6ab1dc0SHong Zhang #undef __FUNCT__
14196da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14206da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1421d6ab1dc0SHong Zhang {
1422d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1423d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1424d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1425d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1426d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1427e745005bSHong Zhang   PetscInt            *adj;
1428e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1429e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1430d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1431ce94432eSBarry Smith   MPI_Comm            comm;
1432d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1433d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1434d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1435d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1436d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1437d6ab1dc0SHong Zhang   MPI_Status          *status;
1438d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1439d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1440a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1441d6ab1dc0SHong Zhang   Mat                 A_loc;
1442d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1443d6ab1dc0SHong Zhang 
1444d6ab1dc0SHong Zhang   PetscFunctionBegin;
1445ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1446d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1447d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1448d6ab1dc0SHong Zhang 
1449d6ab1dc0SHong Zhang   ptap  = c->ptap;
1450d6ab1dc0SHong Zhang   merge = ptap->merge;
1451d6ab1dc0SHong Zhang 
1452e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1453e745005bSHong Zhang   /*------------------------------------------*/
1454d6ab1dc0SHong Zhang   /* get data from symbolic products */
1455d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
14561795a4d1SJed Brown   ierr   = PetscCalloc1((coi[pon]+1),&coa);CHKERRQ(ierr);
1457d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1458d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14591795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1460d6ab1dc0SHong Zhang 
1461d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1462d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1463d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1464d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1465d6ab1dc0SHong Zhang   ai    = a_loc->i;
1466d6ab1dc0SHong Zhang   aj    = a_loc->j;
1467d6ab1dc0SHong Zhang 
1468d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1469d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1470d6ab1dc0SHong Zhang     adj = aj + ai[i];
1471d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1472d6ab1dc0SHong Zhang 
1473d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1474e745005bSHong Zhang     /*-------------------------------------------------------------*/
1475d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1476d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1477d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1478d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1479d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1480d6ab1dc0SHong Zhang       row = poJ[j];
1481d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1482d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1483e745005bSHong Zhang       /* perform sparse axpy */
1484e745005bSHong Zhang       nexta  = 0;
1485d6ab1dc0SHong Zhang       valtmp = pA[j];
1486e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1487e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1488e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1489e745005bSHong Zhang           nexta++;
1490d6ab1dc0SHong Zhang         }
1491e745005bSHong Zhang       }
1492e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1493d6ab1dc0SHong Zhang     }
1494d6ab1dc0SHong Zhang 
1495d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1496d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1497d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1498d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1499d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1500d6ab1dc0SHong Zhang       row = pdJ[j];
1501d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1502d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1503e745005bSHong Zhang       /* perform sparse axpy */
1504e745005bSHong Zhang       nexta  = 0;
1505d6ab1dc0SHong Zhang       valtmp = pA[j];
1506e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1507e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1508e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1509e745005bSHong Zhang           nexta++;
1510d6ab1dc0SHong Zhang         }
1511e745005bSHong Zhang       }
1512e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1513d6ab1dc0SHong Zhang     }
1514d6ab1dc0SHong Zhang   }
1515d6ab1dc0SHong Zhang 
1516d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1517d6ab1dc0SHong Zhang   /*------------------------------------*/
1518d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1519d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1520d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1521d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1522d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1523d6ab1dc0SHong Zhang 
1524dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1525d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1526d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1527d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1528d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1529d6ab1dc0SHong Zhang     k++;
1530d6ab1dc0SHong Zhang   }
1531d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1532d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1533d6ab1dc0SHong Zhang 
1534d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1535d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1536d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1537d6ab1dc0SHong Zhang 
1538d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1539d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1540dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1541d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1542e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1543d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1544d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1545d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1546d6ab1dc0SHong Zhang   }
1547d6ab1dc0SHong Zhang 
1548d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1549d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1550d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1551d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1552d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1553d6ab1dc0SHong Zhang     /* add received vals into ba */
1554d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1555d6ab1dc0SHong Zhang       /* i-th row */
1556d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1557d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1558d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1559d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1560d6ab1dc0SHong Zhang         nextcj = 0;
1561d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1562d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1563d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1564d6ab1dc0SHong Zhang           }
1565d6ab1dc0SHong Zhang         }
1566d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1567d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1568d6ab1dc0SHong Zhang       }
1569d6ab1dc0SHong Zhang     }
1570d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1571d6ab1dc0SHong Zhang   }
1572d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1573d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1574d6ab1dc0SHong Zhang 
1575d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1576d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1577d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1578d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1579d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1580d6ab1dc0SHong Zhang }
1581d6ab1dc0SHong Zhang 
15822bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
15832bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1584d6ab1dc0SHong Zhang #undef __FUNCT__
15856da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
15866da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1587d6ab1dc0SHong Zhang {
1588d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1589d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1590d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15910298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1592d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1593d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1594d6ab1dc0SHong Zhang   PetscInt            nnz;
1595d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1596d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1597ce94432eSBarry Smith   MPI_Comm            comm;
1598d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1599d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1600d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1601d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1602d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1603d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1604d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1605d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1606d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1607d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1608c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1609d6ab1dc0SHong Zhang   PetscScalar         *vals;
1610d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1611d6ab1dc0SHong Zhang 
1612d6ab1dc0SHong Zhang   PetscFunctionBegin;
1613ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1614d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1615d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1616d6ab1dc0SHong 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);
1617d6ab1dc0SHong Zhang   }
1618d6ab1dc0SHong Zhang 
1619d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1620d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1621d6ab1dc0SHong Zhang 
1622d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1623b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1624d6ab1dc0SHong Zhang 
1625d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1626d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16272205254eSKarl Rupp 
1628d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1629d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1630d6ab1dc0SHong Zhang   ai          = a_loc->i;
1631d6ab1dc0SHong Zhang   aj          = a_loc->j;
1632d6ab1dc0SHong Zhang 
1633d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1634d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1635d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1636d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1637d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1638d6ab1dc0SHong Zhang 
1639d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1640d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1641d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1642d6ab1dc0SHong Zhang 
1643d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1644d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1645d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1646785e854fSJed Brown   ierr   = PetscMalloc1((pon+1),&coi);CHKERRQ(ierr);
1647d6ab1dc0SHong Zhang   coi[0] = 0;
1648d6ab1dc0SHong Zhang 
1649d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1650d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1651a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1652d6ab1dc0SHong Zhang   current_space = free_space;
165319f0e57aSHong Zhang 
1654d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1655d6ab1dc0SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1656c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1657d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1658d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1659d6ab1dc0SHong Zhang 
1660d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1661d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1662d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1663d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1664d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1665d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1666d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1667d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1668d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1669d6ab1dc0SHong Zhang     }
1670d6ab1dc0SHong Zhang     nnz = lnk[0];
1671d6ab1dc0SHong Zhang 
1672d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1673d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1674d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1675d6ab1dc0SHong Zhang       nspacedouble++;
1676d6ab1dc0SHong Zhang     }
1677d6ab1dc0SHong Zhang 
1678d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1679d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16802205254eSKarl Rupp 
1681d6ab1dc0SHong Zhang     current_space->array           += nnz;
1682d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1683d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16842205254eSKarl Rupp 
1685d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1686d6ab1dc0SHong Zhang   }
1687d6ab1dc0SHong Zhang 
1688785e854fSJed Brown   ierr = PetscMalloc1((coi[pon]+1),&coj);CHKERRQ(ierr);
1689d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16902205254eSKarl Rupp 
1691118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1692d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1693d6ab1dc0SHong Zhang 
1694d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1695d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1696d6ab1dc0SHong Zhang   /* determine row ownership */
1697b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1698d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
16992205254eSKarl Rupp 
1700d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1701d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17022205254eSKarl Rupp 
1703d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1704d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1705d6ab1dc0SHong Zhang 
1706d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17071795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1708785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17092205254eSKarl Rupp 
1710d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1711d6ab1dc0SHong Zhang   merge->nsend = 0;
1712d6ab1dc0SHong Zhang 
1713785e854fSJed Brown   ierr = PetscMalloc1((size+2),&owners_co);CHKERRQ(ierr);
1714d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1715d6ab1dc0SHong Zhang 
1716d6ab1dc0SHong Zhang   proc = 0;
1717d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1718d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1719d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1720d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1721d6ab1dc0SHong Zhang   }
1722d6ab1dc0SHong Zhang 
1723d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1724d6ab1dc0SHong Zhang   owners_co[0] = 0;
1725d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1726d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1727d6ab1dc0SHong Zhang     if (len_si[proc]) {
1728d6ab1dc0SHong Zhang       merge->nsend++;
1729d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1730d6ab1dc0SHong Zhang       len         += len_si[proc];
1731d6ab1dc0SHong Zhang     }
1732d6ab1dc0SHong Zhang   }
1733d6ab1dc0SHong Zhang 
1734d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17350298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1736d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1737d6ab1dc0SHong Zhang 
1738d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1739d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1740d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1741785e854fSJed Brown   ierr = PetscMalloc1((merge->nsend+1),&swaits);CHKERRQ(ierr);
1742d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1743d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1744d6ab1dc0SHong Zhang     i    = owners_co[proc];
1745d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1746d6ab1dc0SHong Zhang     k++;
1747d6ab1dc0SHong Zhang   }
1748d6ab1dc0SHong Zhang 
1749d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1750785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1751d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1752c280ed6aSJed Brown     PetscMPIInt icompleted;
1753c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1754d6ab1dc0SHong Zhang   }
1755d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1756d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1757d6ab1dc0SHong Zhang 
1758d6ab1dc0SHong Zhang   /* send and recv coi */
1759d6ab1dc0SHong Zhang   /*-------------------*/
1760d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1761d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1762785e854fSJed Brown   ierr   = PetscMalloc1((len+1),&buf_s);CHKERRQ(ierr);
1763d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1764d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1765d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1766d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1767d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1768d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1769d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1770d6ab1dc0SHong Zhang     */
1771d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1772d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1773d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1774d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1775d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1776d6ab1dc0SHong Zhang     nrows       = 0;
1777d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1778d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1779d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1780d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1781d6ab1dc0SHong Zhang       nrows++;
1782d6ab1dc0SHong Zhang     }
1783d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1784d6ab1dc0SHong Zhang     k++;
1785d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1786d6ab1dc0SHong Zhang   }
1787d6ab1dc0SHong Zhang   i = merge->nrecv;
1788d6ab1dc0SHong Zhang   while (i--) {
1789c280ed6aSJed Brown     PetscMPIInt icompleted;
1790c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1791d6ab1dc0SHong Zhang   }
1792d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1793d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1794d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1795d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1796d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1797d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1798d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1799d6ab1dc0SHong Zhang 
1800d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1801d6ab1dc0SHong Zhang   /*------------------------------------------*/
1802d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1803785e854fSJed Brown   ierr  = PetscMalloc1((pn+1),&bi);CHKERRQ(ierr);
1804d6ab1dc0SHong Zhang   bi[0] = 0;
1805d6ab1dc0SHong Zhang 
1806d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1807d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1808a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1809d6ab1dc0SHong Zhang   current_space = free_space;
1810d6ab1dc0SHong Zhang 
1811dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1812d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1813d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1814d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1815d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18162205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1817d6ab1dc0SHong Zhang   }
1818d6ab1dc0SHong Zhang 
1819d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1820d6ab1dc0SHong Zhang   rmax = 0;
1821d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1822d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1823d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1824d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1825d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1826d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1827d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1828d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1829d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1830d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1831d6ab1dc0SHong Zhang     }
1832d6ab1dc0SHong Zhang 
1833d6ab1dc0SHong Zhang     /* add received col data into lnk */
1834d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1835d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1836d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1837d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1838d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1839d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1840d6ab1dc0SHong Zhang       }
1841d6ab1dc0SHong Zhang     }
1842d6ab1dc0SHong Zhang     nnz = lnk[0];
1843d6ab1dc0SHong Zhang 
1844d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1845d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1846d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1847d6ab1dc0SHong Zhang       nspacedouble++;
1848d6ab1dc0SHong Zhang     }
1849d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1850d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1851d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18522205254eSKarl Rupp 
1853d6ab1dc0SHong Zhang     current_space->array           += nnz;
1854d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1855d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18562205254eSKarl Rupp 
1857d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1858d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1859d6ab1dc0SHong Zhang   }
1860d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1861d6ab1dc0SHong Zhang 
1862785e854fSJed Brown   ierr      = PetscMalloc1((bi[pn]+1),&bj);CHKERRQ(ierr);
1863d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1864118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1865d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1866d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1867d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1868d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1869d6ab1dc0SHong Zhang 
1870d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1871d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
18721795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1873d6ab1dc0SHong Zhang 
1874d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1875d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
187633d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1877d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1878d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1879d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1880d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1881d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1882d6ab1dc0SHong Zhang     row  = i + rstart;
1883d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1884d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1885d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1886d6ab1dc0SHong Zhang   }
1887d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1888d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1889d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1890d6ab1dc0SHong Zhang 
1891d6ab1dc0SHong Zhang   merge->bi        = bi;
1892d6ab1dc0SHong Zhang   merge->bj        = bj;
1893d6ab1dc0SHong Zhang   merge->coi       = coi;
1894d6ab1dc0SHong Zhang   merge->coj       = coj;
1895d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1896d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1897d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1898d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1899d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1900d6ab1dc0SHong Zhang 
19016da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1902d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1903d6ab1dc0SHong Zhang 
1904d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1905d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19062205254eSKarl Rupp 
1907d6ab1dc0SHong Zhang   c->ptap     = ptap;
19080298fd71SBarry Smith   ptap->api   = NULL;
19090298fd71SBarry Smith   ptap->apj   = NULL;
1910d6ab1dc0SHong Zhang   ptap->merge = merge;
1911d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
19120298fd71SBarry Smith   ptap->apa   = NULL;
1913f96d37fcSHong Zhang 
1914f96d37fcSHong Zhang   *C = Cmpi;
1915f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1916f96d37fcSHong Zhang   if (bi[pn] != 0) {
191757622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
191857622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1919f96d37fcSHong Zhang   } else {
1920f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1921f96d37fcSHong Zhang   }
1922f96d37fcSHong Zhang #endif
1923f96d37fcSHong Zhang   PetscFunctionReturn(0);
1924f96d37fcSHong Zhang }
1925