xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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;
182499ec78SHong Zhang 
192499ec78SHong Zhang   PetscFunctionBegin;
202499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
21a1a4e84aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
222499ec78SHong Zhang   }
23a2f3521dSMark F. Adams 
24598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
252499ec78SHong Zhang   PetscFunctionReturn(0);
262499ec78SHong Zhang }
272499ec78SHong Zhang 
28f33d1a9aSHong Zhang #undef __FUNCT__
29a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
30a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
31598bc09dSHong Zhang {
32598bc09dSHong Zhang   PetscErrorCode ierr;
33598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
34598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
35598bc09dSHong Zhang 
36598bc09dSHong Zhang   PetscFunctionBegin;
37b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
38598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
39a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
40a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
41ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
42a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
43a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
44598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
45598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
46598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
47598bc09dSHong Zhang   PetscFunctionReturn(0);
48598bc09dSHong Zhang }
49598bc09dSHong Zhang 
502499ec78SHong Zhang #undef __FUNCT__
51a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
52a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
534ae0847bSHong Zhang {
544ae0847bSHong Zhang   PetscErrorCode ierr;
554ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
564ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
574ae0847bSHong Zhang 
584ae0847bSHong Zhang   PetscFunctionBegin;
594ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
6026fbe8dcSKarl Rupp 
614ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
624ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
634ae0847bSHong Zhang   PetscFunctionReturn(0);
644ae0847bSHong Zhang }
654ae0847bSHong Zhang 
664ae0847bSHong Zhang #undef __FUNCT__
67a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
68a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
69598bc09dSHong Zhang {
70598bc09dSHong Zhang   PetscErrorCode ierr;
714ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
72598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
73598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
74598bc09dSHong Zhang   PetscInt       *adi=ad->i,*adj,*aoi=ao->i,*aoj;
75598bc09dSHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
76598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
77598bc09dSHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
78598bc09dSHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
79598bc09dSHong Zhang   PetscInt       cm   =C->rmap->n,anz,pnz;
80598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
8165a57a32SSean Farley   PetscInt       *api,*apj,*apJ,i,j,k,row;
8229825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
83598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
84598bc09dSHong Zhang 
85598bc09dSHong Zhang   PetscFunctionBegin;
86a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
87598bc09dSHong Zhang   /*-----------------------------------------------------*/
88a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
89b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
90a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
91598bc09dSHong Zhang 
92598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
93598bc09dSHong Zhang   /*----------------------------------------------------------*/
94598bc09dSHong Zhang   /* get data from symbolic products */
95a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
96a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
97598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
98598bc09dSHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
99598bc09dSHong Zhang 
100598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
101598bc09dSHong Zhang   apa = ptap->apa;
102598bc09dSHong Zhang 
10329825010SHong Zhang   api = ptap->api;
10429825010SHong Zhang   apj = ptap->apj;
105598bc09dSHong Zhang   for (i=0; i<cm; i++) {
106598bc09dSHong Zhang     /* diagonal portion of A */
107598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
108598bc09dSHong Zhang     adj = ad->j + adi[i];
109598bc09dSHong Zhang     ada = ad->a + adi[i];
110598bc09dSHong Zhang     for (j=0; j<anz; j++) {
111598bc09dSHong Zhang       row = adj[j];
112598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
113598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
114598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
115598bc09dSHong Zhang 
116598bc09dSHong Zhang       /* perform dense axpy */
117598bc09dSHong Zhang       valtmp = ada[j];
118598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
119598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
120598bc09dSHong Zhang       }
121598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
122598bc09dSHong Zhang     }
123598bc09dSHong Zhang 
124598bc09dSHong Zhang     /* off-diagonal portion of A */
125598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
126598bc09dSHong Zhang     aoj = ao->j + aoi[i];
127598bc09dSHong Zhang     aoa = ao->a + aoi[i];
128598bc09dSHong Zhang     for (j=0; j<anz; j++) {
129598bc09dSHong Zhang       row = aoj[j];
130598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
131598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
132598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
133598bc09dSHong Zhang 
134598bc09dSHong Zhang       /* perform dense axpy */
135598bc09dSHong Zhang       valtmp = aoa[j];
136598bc09dSHong Zhang       for (k=0; k<pnz; k++) {
137598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
138598bc09dSHong Zhang       }
139598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
140598bc09dSHong Zhang     }
141598bc09dSHong Zhang 
142598bc09dSHong Zhang     /* set values in C */
143598bc09dSHong Zhang     apJ  = apj + api[i];
144598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
145598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
146598bc09dSHong Zhang 
147598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
148598bc09dSHong Zhang     ca = coa + co->i[i];
149598bc09dSHong Zhang     k  = 0;
150598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
151598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
152598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
153598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
154598bc09dSHong Zhang       k++;
155598bc09dSHong Zhang     }
156598bc09dSHong Zhang 
157598bc09dSHong Zhang     /* diagonal part of C */
158598bc09dSHong Zhang     ca = cda + cd->i[i];
159598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
160598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
161598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
162598bc09dSHong Zhang       k++;
163598bc09dSHong Zhang     }
164598bc09dSHong Zhang 
165598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
166598bc09dSHong Zhang     ca = coa + co->i[i];
167598bc09dSHong Zhang     for (; k0<conz; k0++) {
168598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
169598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
170598bc09dSHong Zhang       k++;
171598bc09dSHong Zhang     }
172598bc09dSHong Zhang   }
173598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
174598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175598bc09dSHong Zhang   PetscFunctionReturn(0);
176598bc09dSHong Zhang }
177598bc09dSHong Zhang 
178598bc09dSHong Zhang #undef __FUNCT__
179a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
180a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
181598bc09dSHong Zhang {
182598bc09dSHong Zhang   PetscErrorCode     ierr;
183*ce94432eSBarry Smith   MPI_Comm           comm;
184bfcd1a73SHong Zhang   Mat                Cmpi;
185598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1860298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1874ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
188bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1894ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
1904ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
191d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
192bfcd1a73SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
193598bc09dSHong Zhang   PetscBT            lnkbt;
194598bc09dSHong Zhang   PetscScalar        *apa;
195bfcd1a73SHong Zhang   PetscReal          afill;
196c397d730SHong Zhang   PetscBool          scalable=PETSC_TRUE;
197f84c536bSHong Zhang   PetscInt           nlnk_max,armax,prmax;
198598bc09dSHong Zhang 
199598bc09dSHong Zhang   PetscFunctionBegin;
200*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
201a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
202cf1a0672SHong 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);
203a1a4e84aSHong Zhang   }
204152983bfSHong Zhang 
205152983bfSHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
2060298fd71SBarry Smith   ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,NULL);CHKERRQ(ierr);
207cf1a0672SHong Zhang   if (scalable) {
2085c913ed7SHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);CHKERRQ(ierr);
2091634b032SHong Zhang     PetscFunctionReturn(0);
2101634b032SHong Zhang   }
211152983bfSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
212a1a4e84aSHong Zhang 
213598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
214598bc09dSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
215598bc09dSHong Zhang 
216598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
217b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
21819f0e57aSHong Zhang 
219598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
220a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
221598bc09dSHong Zhang 
222a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
223a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
224598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
225598bc09dSHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
226598bc09dSHong Zhang 
227598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
228598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
229598bc09dSHong Zhang   ierr      = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
230a1a4e84aSHong Zhang   ptap->api = api;
231598bc09dSHong Zhang   api[0]    = 0;
232598bc09dSHong Zhang 
233598bc09dSHong Zhang   /* create and initialize a linked list */
234f84c536bSHong Zhang   armax    = ad->rmax+ao->rmax;
235f84c536bSHong Zhang   prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
236f84c536bSHong Zhang   nlnk_max = armax*prmax;
237f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2380d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
239598bc09dSHong Zhang 
240bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
241bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
2422205254eSKarl Rupp 
243598bc09dSHong Zhang   current_space = free_space;
244598bc09dSHong Zhang 
245598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
246598bc09dSHong Zhang   for (i=0; i<am; i++) {
247598bc09dSHong Zhang     /* diagonal portion of A */
248598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
249598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
250598bc09dSHong Zhang       row  = *adj++;
251598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
252598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
253598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2541fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
255598bc09dSHong Zhang     }
256598bc09dSHong Zhang     /* off-diagonal portion of A */
257598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
258598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
259598bc09dSHong Zhang       row  = *aoj++;
260598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
261598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2621fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
263598bc09dSHong Zhang     }
264598bc09dSHong Zhang 
265d14fa243SHong Zhang     apnz     = lnk[0];
266598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
267598bc09dSHong Zhang 
268598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
269598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
270598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
271598bc09dSHong Zhang       nspacedouble++;
272598bc09dSHong Zhang     }
273598bc09dSHong Zhang 
274598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2751fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
276598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2772205254eSKarl Rupp 
278598bc09dSHong Zhang     current_space->array           += apnz;
279598bc09dSHong Zhang     current_space->local_used      += apnz;
280598bc09dSHong Zhang     current_space->local_remaining -= apnz;
281598bc09dSHong Zhang   }
282598bc09dSHong Zhang 
283598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
284598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
285a1a4e84aSHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
286a1a4e84aSHong Zhang   apj  = ptap->apj;
287a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
288598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
289598bc09dSHong Zhang 
290f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
291f84c536bSHong Zhang   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
292f84c536bSHong Zhang   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
2932205254eSKarl Rupp 
294f84c536bSHong Zhang   ptap->apa = apa;
295f84c536bSHong Zhang 
296bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
297598bc09dSHong Zhang   /*----------------------------------------------------*/
298bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
299bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
300a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,A->rmap->bs,P->cmap->bs);CHKERRQ(ierr);
301a2f3521dSMark F. Adams 
302bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
303bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
304598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
305598bc09dSHong Zhang   for (i=0; i<am; i++) {
306598bc09dSHong Zhang     row  = i + rstart;
307598bc09dSHong Zhang     apnz = api[i+1] - api[i];
308bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
309598bc09dSHong Zhang     apj += apnz;
310598bc09dSHong Zhang   }
311bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
312bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313598bc09dSHong Zhang 
314bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
315bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
316bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
317bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
318598bc09dSHong Zhang 
319bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
320bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
321598bc09dSHong Zhang   c->ptap = ptap;
322598bc09dSHong Zhang 
323bfcd1a73SHong Zhang   *C = Cmpi;
324bfcd1a73SHong Zhang 
325bfcd1a73SHong Zhang   /* set MatInfo */
326118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
327bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
328bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
329bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
330bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
331bfcd1a73SHong Zhang 
332bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
333bfcd1a73SHong Zhang   if (api[am]) {
334bfcd1a73SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
335bfcd1a73SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
336bfcd1a73SHong Zhang   } else {
337bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
338bfcd1a73SHong Zhang   }
339bfcd1a73SHong Zhang #endif
340598bc09dSHong Zhang   PetscFunctionReturn(0);
341598bc09dSHong Zhang }
342598bc09dSHong Zhang 
3439123193aSHong Zhang #undef __FUNCT__
3449123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3459123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3469123193aSHong Zhang {
3479123193aSHong Zhang   PetscErrorCode ierr;
3489123193aSHong Zhang 
3499123193aSHong Zhang   PetscFunctionBegin;
3509123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3519123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3529123193aSHong Zhang   }
3539123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3549123193aSHong Zhang   PetscFunctionReturn(0);
3559123193aSHong Zhang }
3569123193aSHong Zhang 
3574324174eSBarry Smith typedef struct {
3584324174eSBarry Smith   Mat         workB;
3594324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3604324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3614324174eSBarry Smith } MPIAIJ_MPIDense;
3624324174eSBarry Smith 
3637af9e4fdSHong Zhang #undef __FUNCT__
36496b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
36596b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3664324174eSBarry Smith {
3674324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3684324174eSBarry Smith   PetscErrorCode  ierr;
3694324174eSBarry Smith 
3704324174eSBarry Smith   PetscFunctionBegin;
3716bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3724324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3734324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3744324174eSBarry Smith   PetscFunctionReturn(0);
3754324174eSBarry Smith }
3764324174eSBarry Smith 
3779123193aSHong Zhang #undef __FUNCT__
3789123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
3799123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
3809123193aSHong Zhang {
3819123193aSHong Zhang   PetscErrorCode         ierr;
382f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
383d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
384bf0cc555SLisandro Dalcin   PetscContainer         container;
3854324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
3864324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
3874324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
3884324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
389d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
3909123193aSHong Zhang 
3919123193aSHong Zhang   PetscFunctionBegin;
392*ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
393d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
394a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(*C,A->rmap->bs,B->cmap->bs);CHKERRQ(ierr);
395cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
3960298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
397cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
398cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3992205254eSKarl Rupp 
400d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
401f32d5d43SBarry Smith 
4024324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
403f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4040298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4054324174eSBarry Smith   /* Create work arrays needed */
406d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
407d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
4084324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
4094324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
4104324174eSBarry Smith 
411*ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
412bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
41396b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
414bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
415bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4169123193aSHong Zhang   PetscFunctionReturn(0);
4179123193aSHong Zhang }
4189123193aSHong Zhang 
4197af9e4fdSHong Zhang #undef __FUNCT__
4207af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
421f32d5d43SBarry Smith /*
4222636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4232636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
424f32d5d43SBarry Smith */
4254324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
426f32d5d43SBarry Smith {
427f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
428f32d5d43SBarry Smith   PetscErrorCode         ierr;
429f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
430f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
431f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
432f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
433f32d5d43SBarry Smith   PetscInt               i,j,k;
434aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
435aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
436f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
437*ce94432eSBarry Smith   MPI_Comm               comm;
438d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
439f32d5d43SBarry Smith   MPI_Status             status;
4404324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
441bf0cc555SLisandro Dalcin   PetscContainer         container;
4424324174eSBarry Smith   Mat                    workB;
443f32d5d43SBarry Smith 
444f32d5d43SBarry Smith   PetscFunctionBegin;
445*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
446bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
44729825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
448bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4494324174eSBarry Smith 
4504324174eSBarry Smith   workB = *outworkB = contents->workB;
451cf1a0672SHong 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);
452f32d5d43SBarry Smith   sindices = to->indices;
453f32d5d43SBarry Smith   sstarts  = to->starts;
454f32d5d43SBarry Smith   sprocs   = to->procs;
4554324174eSBarry Smith   swaits   = contents->swaits;
4564324174eSBarry Smith   svalues  = contents->svalues;
457f32d5d43SBarry Smith 
458f32d5d43SBarry Smith   rindices = from->indices;
459f32d5d43SBarry Smith   rstarts  = from->starts;
460f32d5d43SBarry Smith   rprocs   = from->procs;
4614324174eSBarry Smith   rwaits   = contents->rwaits;
4624324174eSBarry Smith   rvalues  = contents->rvalues;
463f32d5d43SBarry Smith 
4648c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
4658c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
466f32d5d43SBarry Smith 
467f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
468f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
469f32d5d43SBarry Smith   }
470f32d5d43SBarry Smith 
471f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
472f32d5d43SBarry Smith     /* pack a message at a time */
473f32d5d43SBarry Smith     CHKMEMQ;
474f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
475f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
476f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
4772636ff26SBarry Smith       }
4782636ff26SBarry Smith     }
479f32d5d43SBarry Smith     CHKMEMQ;
480f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
481f32d5d43SBarry Smith   }
4822636ff26SBarry Smith 
483f32d5d43SBarry Smith   nrecvs = from->n;
484f32d5d43SBarry Smith   while (nrecvs) {
485f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
486f32d5d43SBarry Smith     nrecvs--;
487f32d5d43SBarry Smith     /* unpack a message at a time */
488f32d5d43SBarry Smith     CHKMEMQ;
489f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
490f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
491f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
4922636ff26SBarry Smith       }
4932636ff26SBarry Smith     }
494f32d5d43SBarry Smith     CHKMEMQ;
495f32d5d43SBarry Smith   }
496cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
497f32d5d43SBarry Smith 
4988c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
4998c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
500f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
501f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
502f32d5d43SBarry Smith   PetscFunctionReturn(0);
503f32d5d43SBarry Smith }
504f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
505f32d5d43SBarry Smith 
5069123193aSHong Zhang #undef __FUNCT__
5079123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5089123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5099123193aSHong Zhang {
5109123193aSHong Zhang   PetscErrorCode ierr;
511f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
512f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
513f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
514f32d5d43SBarry Smith   Mat            workB;
5159123193aSHong Zhang 
5169123193aSHong Zhang   PetscFunctionBegin;
517f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
518f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
519f32d5d43SBarry Smith 
520f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5214324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
522f32d5d43SBarry Smith 
523f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
524f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5259123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5269123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5279123193aSHong Zhang   PetscFunctionReturn(0);
5289123193aSHong Zhang }
529cf1a0672SHong Zhang 
5301634b032SHong Zhang #undef __FUNCT__
531cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
532cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
5331634b032SHong Zhang {
5341634b032SHong Zhang   PetscErrorCode ierr;
535cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
536cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
537cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
538cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
539cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
540cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
541cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
54229825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
54329825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
544cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
54529825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
546cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
54729825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
54829825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
5491634b032SHong Zhang 
5501634b032SHong Zhang   PetscFunctionBegin;
551cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
552cf1a0672SHong Zhang   /*-----------------------------------------------------*/
553cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
554b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
555cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
556cf1a0672SHong Zhang 
557cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
558cf1a0672SHong Zhang   /*----------------------------------------------------------*/
559cf1a0672SHong Zhang   /* get data from symbolic products */
560cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
561cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
562cf1a0672SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
563cf1a0672SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
564cf1a0672SHong Zhang 
56529825010SHong Zhang   api = ptap->api;
56629825010SHong Zhang   apj = ptap->apj;
567cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
56829825010SHong Zhang     apJ = apj + api[i];
56929825010SHong Zhang 
570cf1a0672SHong Zhang     /* diagonal portion of A */
571cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
572cf1a0672SHong Zhang     adj = ad->j + adi[i];
573cf1a0672SHong Zhang     ada = ad->a + adi[i];
574cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
575cf1a0672SHong Zhang       row = adj[j];
576cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
577cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
578cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
57929825010SHong Zhang       /* perform sparse axpy */
580cf1a0672SHong Zhang       valtmp = ada[j];
58129825010SHong Zhang       nextp  = 0;
58229825010SHong Zhang       for (k=0; nextp<pnz; k++) {
58329825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
58429825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
58529825010SHong Zhang         }
5861634b032SHong Zhang       }
587cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
588cf1a0672SHong Zhang     }
5891634b032SHong Zhang 
590cf1a0672SHong Zhang     /* off-diagonal portion of A */
591cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
592cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
593cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
594cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
595cf1a0672SHong Zhang       row = aoj[j];
596cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
597cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
598cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
59929825010SHong Zhang       /* perform sparse axpy */
600cf1a0672SHong Zhang       valtmp = aoa[j];
60129825010SHong Zhang       nextp  = 0;
60229825010SHong Zhang       for (k=0; nextp<pnz; k++) {
60329825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
60429825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
60529825010SHong Zhang         }
606cf1a0672SHong Zhang       }
607cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
608cf1a0672SHong Zhang     }
6091634b032SHong Zhang 
610cf1a0672SHong Zhang     /* set values in C */
611cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
612cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6131634b032SHong Zhang 
614cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
615cf1a0672SHong Zhang     ca = coa + co->i[i];
616cf1a0672SHong Zhang     k  = 0;
617cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
618cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
61929825010SHong Zhang       ca[k0]        = apa_sparse[k];
62029825010SHong Zhang       apa_sparse[k] = 0.0;
621cf1a0672SHong Zhang       k++;
622cf1a0672SHong Zhang     }
6231634b032SHong Zhang 
624cf1a0672SHong Zhang     /* diagonal part of C */
625cf1a0672SHong Zhang     ca = cda + cd->i[i];
626cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
62729825010SHong Zhang       ca[k1]        = apa_sparse[k];
62829825010SHong Zhang       apa_sparse[k] = 0.0;
629cf1a0672SHong Zhang       k++;
630cf1a0672SHong Zhang     }
631cf1a0672SHong Zhang 
632cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
633cf1a0672SHong Zhang     ca = coa + co->i[i];
634cf1a0672SHong Zhang     for (; k0<conz; k0++) {
63529825010SHong Zhang       ca[k0]        = apa_sparse[k];
63629825010SHong Zhang       apa_sparse[k] = 0.0;
637cf1a0672SHong Zhang       k++;
638cf1a0672SHong Zhang     }
639cf1a0672SHong Zhang   }
640cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
641cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
642cf1a0672SHong Zhang   PetscFunctionReturn(0);
643cf1a0672SHong Zhang }
644cf1a0672SHong Zhang 
645cf1a0672SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */
646cf1a0672SHong Zhang #undef __FUNCT__
647cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
648cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C)
649cf1a0672SHong Zhang {
650cf1a0672SHong Zhang   PetscErrorCode     ierr;
651*ce94432eSBarry Smith   MPI_Comm           comm;
652cf1a0672SHong Zhang   Mat                Cmpi;
653cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
6540298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
655cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
656cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
657cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
658cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
659f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
660cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
661cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
662cf1a0672SHong Zhang   PetscReal          afill;
663cf1a0672SHong Zhang   PetscScalar        *apa;
664cf1a0672SHong Zhang 
665cf1a0672SHong Zhang   PetscFunctionBegin;
666*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
667cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
668cf1a0672SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
669cf1a0672SHong Zhang 
670cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
671b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
67219f0e57aSHong Zhang 
673cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
674cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
675cf1a0672SHong Zhang 
676cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
677cf1a0672SHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
678cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
679cf1a0672SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
680cf1a0672SHong Zhang 
681cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
682cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
683cf1a0672SHong Zhang   ierr      = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
684cf1a0672SHong Zhang   ptap->api = api;
685cf1a0672SHong Zhang   api[0]    = 0;
686cf1a0672SHong Zhang 
687cf1a0672SHong Zhang   /* create and initialize a linked list */
688cf1a0672SHong Zhang   armax    = ad->rmax+ao->rmax;
689cf1a0672SHong Zhang   prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
690cf1a0672SHong Zhang   nlnk_max = armax*prmax;
691cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
692f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
693cf1a0672SHong Zhang 
694cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
695cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
6962205254eSKarl Rupp 
697cf1a0672SHong Zhang   current_space = free_space;
698cf1a0672SHong Zhang 
699cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
700cf1a0672SHong Zhang   for (i=0; i<am; i++) {
701cf1a0672SHong Zhang     /* diagonal portion of A */
702cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
703cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
704cf1a0672SHong Zhang       row  = *adj++;
705cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
706cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
707cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
708f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
709cf1a0672SHong Zhang     }
710cf1a0672SHong Zhang     /* off-diagonal portion of A */
711cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
712cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
713cf1a0672SHong Zhang       row  = *aoj++;
714cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
715cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
716f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
717cf1a0672SHong Zhang     }
718cf1a0672SHong Zhang 
719f84c536bSHong Zhang     apnz     = *lnk;
720cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
721cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
722cf1a0672SHong Zhang 
723cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
724cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
725cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
726cf1a0672SHong Zhang       nspacedouble++;
727cf1a0672SHong Zhang     }
728cf1a0672SHong Zhang 
729cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
730f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
731cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
7322205254eSKarl Rupp 
733cf1a0672SHong Zhang     current_space->array           += apnz;
734cf1a0672SHong Zhang     current_space->local_used      += apnz;
735cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
736cf1a0672SHong Zhang   }
737cf1a0672SHong Zhang 
738cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
739cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
740cf1a0672SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
741cf1a0672SHong Zhang   apj  = ptap->apj;
742cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
743f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
744cf1a0672SHong Zhang 
745cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
746cf1a0672SHong Zhang   /*----------------------------------------------------*/
747cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
748cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
749a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,A->rmap->bs,P->cmap->bs);CHKERRQ(ierr);
750cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
751cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
752cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
753cf1a0672SHong Zhang 
75429825010SHong Zhang   /* malloc apa for assembly Cmpi */
755cf1a0672SHong Zhang   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
756cf1a0672SHong Zhang   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
7572205254eSKarl Rupp 
75829825010SHong Zhang   ptap->apa = apa;
759cf1a0672SHong Zhang   for (i=0; i<am; i++) {
760cf1a0672SHong Zhang     row  = i + rstart;
761cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
762cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
763cf1a0672SHong Zhang     apj += apnz;
764cf1a0672SHong Zhang   }
765cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
766cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767cf1a0672SHong Zhang 
768cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
769cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
770cf1a0672SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
771cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
772cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
773cf1a0672SHong Zhang 
774cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
775cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
776cf1a0672SHong Zhang   c->ptap = ptap;
777cf1a0672SHong Zhang 
778cf1a0672SHong Zhang   *C = Cmpi;
779cf1a0672SHong Zhang 
780cf1a0672SHong Zhang   /* set MatInfo */
781118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
782cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
783cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
784cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
785cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
786cf1a0672SHong Zhang 
787cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
788cf1a0672SHong Zhang   if (api[am]) {
789cf1a0672SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
790cf1a0672SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
791cf1a0672SHong Zhang   } else {
792cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
793cf1a0672SHong Zhang   }
794cf1a0672SHong Zhang #endif
7951634b032SHong Zhang   PetscFunctionReturn(0);
7961634b032SHong Zhang }
797f96d37fcSHong Zhang 
798f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
799f96d37fcSHong Zhang #undef __FUNCT__
800f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
801c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
802f96d37fcSHong Zhang {
803f96d37fcSHong Zhang   PetscErrorCode ierr;
804c397d730SHong Zhang   PetscBool      scalable=PETSC_TRUE,viamatmatmult=PETSC_FALSE;
805f96d37fcSHong Zhang 
806f96d37fcSHong Zhang   PetscFunctionBegin;
8070298fd71SBarry Smith   ierr = PetscOptionsBool("-mattransposematmult_viamatmatmult","Use R=Pt and C=R*A","",viamatmatmult,&viamatmatmult,NULL);CHKERRQ(ierr);
808ab1b013aSHong Zhang   if (viamatmatmult) {
809ab1b013aSHong Zhang     Mat         Pt;
810ab1b013aSHong Zhang     Mat_PtAPMPI *ptap;
811ab1b013aSHong Zhang     Mat_MPIAIJ  *c;
812ab1b013aSHong Zhang     if (scall == MAT_INITIAL_MATRIX) {
813ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
814ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
815ab1b013aSHong Zhang 
816ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
817ab1b013aSHong Zhang       ptap     = c->ptap;
818ab1b013aSHong Zhang       ptap->Pt = Pt;
819ab1b013aSHong Zhang     } else if (scall == MAT_REUSE_MATRIX) {
820ab1b013aSHong Zhang       c    = (Mat_MPIAIJ*)(*C)->data;
821ab1b013aSHong Zhang       ptap = c->ptap;
822ab1b013aSHong Zhang       Pt   = ptap->Pt;
823ab1b013aSHong Zhang       ierr = MatTranspose(P,scall,&Pt);CHKERRQ(ierr);
824ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,scall,fill,C);CHKERRQ(ierr);
825*ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not supported");
826ab1b013aSHong Zhang     PetscFunctionReturn(0);
827ab1b013aSHong Zhang   }
828ab1b013aSHong Zhang 
829f96d37fcSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
830d6ab1dc0SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
8310298fd71SBarry Smith     ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,NULL);CHKERRQ(ierr);
832d6ab1dc0SHong Zhang     if  (scalable) {
833d6ab1dc0SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr);
834d6ab1dc0SHong Zhang     } else {
835c5af039cSHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
836f96d37fcSHong Zhang     }
837d6ab1dc0SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
838d6ab1dc0SHong Zhang   }
839d6ab1dc0SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
840f96d37fcSHong Zhang   PetscFunctionReturn(0);
841f96d37fcSHong Zhang }
842f96d37fcSHong Zhang 
843f96d37fcSHong Zhang #undef __FUNCT__
844f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
845c5af039cSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
846f96d37fcSHong Zhang {
847c5af039cSHong Zhang   PetscErrorCode      ierr;
848c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
849497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
850c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
851c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
852d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
853497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
854e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
855c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
856*ce94432eSBarry Smith   MPI_Comm            comm;
857c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
858c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
859c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
860c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
861c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
862c5af039cSHong Zhang   MPI_Status          *status;
863c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
864d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
865a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
866c5af039cSHong Zhang   Mat                 A_loc;
867c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
868f96d37fcSHong Zhang 
869f96d37fcSHong Zhang   PetscFunctionBegin;
870*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
871c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
872c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
873c5af039cSHong Zhang 
874c5af039cSHong Zhang   ptap  = c->ptap;
875c5af039cSHong Zhang   merge = ptap->merge;
876c5af039cSHong Zhang 
877c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
878c5af039cSHong Zhang   /*--------------------------------------------------------------*/
879c5af039cSHong Zhang   /* get data from symbolic products */
880c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
881c5af039cSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
882c5af039cSHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
883c5af039cSHong Zhang 
884c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
885c5af039cSHong Zhang   owners = merge->rowmap->range;
886c5af039cSHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
887c5af039cSHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
888c5af039cSHong Zhang 
889c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
890c5af039cSHong Zhang   A_loc = ptap->A_loc;
891c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
892c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
893d6ab1dc0SHong Zhang   ai    = a_loc->i;
894d6ab1dc0SHong Zhang   aj    = a_loc->j;
895c5af039cSHong Zhang 
896e745005bSHong Zhang   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */
897e745005bSHong Zhang   ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
898c5af039cSHong Zhang 
899c5af039cSHong Zhang   for (i=0; i<am; i++) {
900e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
901d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
902d6ab1dc0SHong Zhang     adj = aj + ai[i];
903d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
904c5af039cSHong Zhang     for (j=0; j<anz; j++) {
905e745005bSHong Zhang       aval[adj[j]] = ada[j];
906c5af039cSHong Zhang     }
907c5af039cSHong Zhang 
908c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
909c5af039cSHong Zhang     /*--------------------------------------------------------------*/
910c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
911c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
912c5af039cSHong Zhang     poJ = po->j + po->i[i];
913c5af039cSHong Zhang     pA  = po->a + po->i[i];
914c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
915c5af039cSHong Zhang       row = poJ[j];
916c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
917c5af039cSHong Zhang       cj  = coj + coi[row];
918c5af039cSHong Zhang       ca  = coa + coi[row];
919c5af039cSHong Zhang       /* perform dense axpy */
920c5af039cSHong Zhang       valtmp = pA[j];
921c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
922e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
923c5af039cSHong Zhang       }
924c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
925c5af039cSHong Zhang     }
926c5af039cSHong Zhang 
927c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
928c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
929c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
930c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
931c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
932c5af039cSHong Zhang       row = pdJ[j];
933c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
934c5af039cSHong Zhang       cj  = bj + bi[row];
935c5af039cSHong Zhang       ca  = ba + bi[row];
936c5af039cSHong Zhang       /* perform dense axpy */
937c5af039cSHong Zhang       valtmp = pA[j];
938c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
939e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
940c5af039cSHong Zhang       }
941c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
942c5af039cSHong Zhang     }
943c5af039cSHong Zhang 
944d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
945d6ab1dc0SHong Zhang     aJ = aj + ai[i];
946e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
947c5af039cSHong Zhang   }
948c5af039cSHong Zhang 
949c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
950c5af039cSHong Zhang   /*------------------------------------*/
951c5af039cSHong Zhang   buf_ri = merge->buf_ri;
952c5af039cSHong Zhang   buf_rj = merge->buf_rj;
953c5af039cSHong Zhang   len_s  = merge->len_s;
954c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
955c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
956c5af039cSHong Zhang 
957c5af039cSHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
958c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
959c5af039cSHong Zhang     if (!len_s[proc]) continue;
960c5af039cSHong Zhang     i    = merge->owners_co[proc];
961c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
962c5af039cSHong Zhang     k++;
963c5af039cSHong Zhang   }
964c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
965c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
966c5af039cSHong Zhang 
967c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
968c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
969c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
970c5af039cSHong Zhang 
971c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
972c5af039cSHong Zhang   /*----------------------------------------------------*/
973c5af039cSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
974c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
975c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
976c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
977c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
978c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
979c5af039cSHong Zhang   }
980c5af039cSHong Zhang 
981c5af039cSHong Zhang   for (i=0; i<cm; i++) {
982c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
983c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
984c5af039cSHong Zhang     ba_i = ba + bi[i];
985c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
986c5af039cSHong Zhang     /* add received vals into ba */
987c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
988c5af039cSHong Zhang       /* i-th row */
989c5af039cSHong Zhang       if (i == *nextrow[k]) {
990c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
991c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
992c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
993c5af039cSHong Zhang         nextcj = 0;
994c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
995c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
996c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
997c5af039cSHong Zhang           }
998c5af039cSHong Zhang         }
999c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1000c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1001c5af039cSHong Zhang       }
1002c5af039cSHong Zhang     }
1003c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1004c5af039cSHong Zhang   }
1005c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1006c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1007c5af039cSHong Zhang 
1008c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1009c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1010c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1011c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1012e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1013f96d37fcSHong Zhang   PetscFunctionReturn(0);
1014f96d37fcSHong Zhang }
1015f96d37fcSHong Zhang 
1016c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1017f96d37fcSHong Zhang #undef __FUNCT__
1018f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
1019f96d37fcSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1020f96d37fcSHong Zhang {
1021f96d37fcSHong Zhang   PetscErrorCode      ierr;
10224e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1023f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
10240298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1025f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1026f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1027f96d37fcSHong Zhang   PetscInt            nnz;
10284e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1029497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1030f96d37fcSHong Zhang   PetscBT             lnkbt;
1031*ce94432eSBarry Smith   MPI_Comm            comm;
1032f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1033f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1034f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1035f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1036f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1037f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1038f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1039f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1040d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1041f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1042c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1043f96d37fcSHong Zhang   PetscScalar         *vals;
10444e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1045f96d37fcSHong Zhang 
1046f96d37fcSHong Zhang   PetscFunctionBegin;
1047*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1048c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1049c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1050c5af039cSHong 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);
1051c5af039cSHong Zhang   }
1052c5af039cSHong Zhang 
1053f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1054f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1055f96d37fcSHong Zhang 
1056f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1057f96d37fcSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1058f96d37fcSHong Zhang 
1059f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1060f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
10612205254eSKarl Rupp 
1062c5af039cSHong Zhang   ptap->A_loc = A_loc;
10632205254eSKarl Rupp 
10641c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1065d6ab1dc0SHong Zhang   ai    = a_loc->i;
1066d6ab1dc0SHong Zhang   aj    = a_loc->j;
1067f96d37fcSHong Zhang 
1068f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1069f96d37fcSHong Zhang   /*----------------------------------------------------*/
10704e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
10714e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
10724e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
10734e938277SHong Zhang 
10744e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
10754e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
10764e938277SHong Zhang   poti = pot->i; potj = pot->j;
1077f96d37fcSHong Zhang 
1078f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
10792205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1080f96d37fcSHong Zhang   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1081f96d37fcSHong Zhang   coi[0] = 0;
1082f96d37fcSHong Zhang 
1083f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1084d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1085a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1086f96d37fcSHong Zhang   current_space = free_space;
108719f0e57aSHong Zhang 
1088f96d37fcSHong Zhang   /* create and initialize a linked list */
10894e938277SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1090c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
10914e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
10924e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1093f96d37fcSHong Zhang 
1094f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1095f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1096f96d37fcSHong Zhang     ptJ = potj + poti[i];
1097f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1098f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1099d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1100d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1101f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1102d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1103f96d37fcSHong Zhang     }
11044e938277SHong Zhang     nnz = lnk[0];
1105f96d37fcSHong Zhang 
1106f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1107f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1108f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1109f96d37fcSHong Zhang       nspacedouble++;
1110f96d37fcSHong Zhang     }
1111f96d37fcSHong Zhang 
1112f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11134e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11142205254eSKarl Rupp 
1115f96d37fcSHong Zhang     current_space->array           += nnz;
1116f96d37fcSHong Zhang     current_space->local_used      += nnz;
1117f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11182205254eSKarl Rupp 
1119f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1120f96d37fcSHong Zhang   }
1121f96d37fcSHong Zhang 
1122f96d37fcSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1123f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
11242205254eSKarl Rupp 
1125118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1126f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1127f96d37fcSHong Zhang 
1128f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1129f96d37fcSHong Zhang   /*----------------------------------------------*/
1130f96d37fcSHong Zhang   /* determine row ownership */
1131f96d37fcSHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1132f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11332205254eSKarl Rupp 
1134f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1135f96d37fcSHong Zhang   merge->rowmap->bs = 1;
11362205254eSKarl Rupp 
1137f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1138f96d37fcSHong Zhang   owners = merge->rowmap->range;
1139f96d37fcSHong Zhang 
1140f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
1141f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1142f96d37fcSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1143f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
11442205254eSKarl Rupp 
1145f96d37fcSHong Zhang   len_s        = merge->len_s;
1146f96d37fcSHong Zhang   merge->nsend = 0;
1147f96d37fcSHong Zhang 
1148f96d37fcSHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1149f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1150f96d37fcSHong Zhang 
1151f96d37fcSHong Zhang   proc = 0;
1152f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1153f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1154f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1155f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1156f96d37fcSHong Zhang   }
1157f96d37fcSHong Zhang 
1158f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1159f96d37fcSHong Zhang   owners_co[0] = 0;
1160f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1161f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1162f96d37fcSHong Zhang     if (len_si[proc]) {
1163f96d37fcSHong Zhang       merge->nsend++;
1164f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1165f96d37fcSHong Zhang       len         += len_si[proc];
1166f96d37fcSHong Zhang     }
1167f96d37fcSHong Zhang   }
1168f96d37fcSHong Zhang 
1169f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
11700298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1171f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1172f96d37fcSHong Zhang 
1173f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1174f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1175f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1176f96d37fcSHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1177f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1178f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1179f96d37fcSHong Zhang     i    = owners_co[proc];
1180f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1181f96d37fcSHong Zhang     k++;
1182f96d37fcSHong Zhang   }
1183f96d37fcSHong Zhang 
1184f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1185f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1186f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1187c280ed6aSJed Brown     PetscMPIInt icompleted;
1188c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1189f96d37fcSHong Zhang   }
1190f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1191f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1192f96d37fcSHong Zhang 
1193f96d37fcSHong Zhang   /* send and recv coi */
1194f96d37fcSHong Zhang   /*-------------------*/
1195f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1196f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1197f96d37fcSHong Zhang   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1198f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1199f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1200f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1201f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1202f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1203f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1204f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1205f96d37fcSHong Zhang     */
1206f96d37fcSHong Zhang     /*-------------------------------------------*/
1207f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1208f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1209f96d37fcSHong Zhang     buf_si[0]   = nrows;
1210f96d37fcSHong Zhang     buf_si_i[0] = 0;
1211f96d37fcSHong Zhang     nrows       = 0;
1212f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1213f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1214f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1215f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1216f96d37fcSHong Zhang       nrows++;
1217f96d37fcSHong Zhang     }
1218f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1219f96d37fcSHong Zhang     k++;
1220f96d37fcSHong Zhang     buf_si += len_si[proc];
1221f96d37fcSHong Zhang   }
1222f96d37fcSHong Zhang   i = merge->nrecv;
1223f96d37fcSHong Zhang   while (i--) {
1224c280ed6aSJed Brown     PetscMPIInt icompleted;
1225c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1226f96d37fcSHong Zhang   }
1227f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1228f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1229f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1230f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1231f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1232f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1233f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1234f96d37fcSHong Zhang 
1235f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1236f96d37fcSHong Zhang   /*------------------------------------------*/
1237f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1238f96d37fcSHong Zhang   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1239f96d37fcSHong Zhang   bi[0] = 0;
1240f96d37fcSHong Zhang 
1241c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1242d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1243a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1244f96d37fcSHong Zhang   current_space = free_space;
1245f96d37fcSHong Zhang 
1246f96d37fcSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1247f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1248f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1249f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1250f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1251f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1252f96d37fcSHong Zhang   }
1253f96d37fcSHong Zhang 
12541c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1255f96d37fcSHong Zhang   rmax = 0;
1256f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1257f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1258f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1259f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1260f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1261f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1262d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1263d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1264f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1265d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1266f96d37fcSHong Zhang     }
12674e938277SHong Zhang 
1268f96d37fcSHong Zhang     /* add received col data into lnk */
1269f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1270f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1271f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1272f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
12734e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1274f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1275f96d37fcSHong Zhang       }
1276f96d37fcSHong Zhang     }
12774e938277SHong Zhang     nnz = lnk[0];
1278f96d37fcSHong Zhang 
1279f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1280f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1281f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1282f96d37fcSHong Zhang       nspacedouble++;
1283f96d37fcSHong Zhang     }
1284f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
12854e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1286f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
12872205254eSKarl Rupp 
1288f96d37fcSHong Zhang     current_space->array           += nnz;
1289f96d37fcSHong Zhang     current_space->local_used      += nnz;
1290f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
12912205254eSKarl Rupp 
1292f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1293f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1294f96d37fcSHong Zhang   }
1295f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1296f96d37fcSHong Zhang 
1297f96d37fcSHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1298f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
12992205254eSKarl Rupp 
1300118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1301f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1302d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13034e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13044e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1305f96d37fcSHong Zhang 
13061c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13071c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13081c7d5954SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
13091c7d5954SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1310f96d37fcSHong Zhang 
1311f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1312f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1313a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1314f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1315f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1316f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1317f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1318f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1319f96d37fcSHong Zhang     row  = i + rstart;
1320f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1321f96d37fcSHong Zhang     Jptr = bj + bi[i];
1322f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1323f96d37fcSHong Zhang   }
1324f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1325f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1326f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1327f96d37fcSHong Zhang 
1328f96d37fcSHong Zhang   merge->bi        = bi;
1329f96d37fcSHong Zhang   merge->bj        = bj;
1330f96d37fcSHong Zhang   merge->coi       = coi;
1331f96d37fcSHong Zhang   merge->coj       = coj;
1332f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1333f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1334f96d37fcSHong Zhang   merge->owners_co = owners_co;
1335f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1336f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1337f96d37fcSHong Zhang 
1338d6ab1dc0SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1339c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1340f96d37fcSHong Zhang 
1341f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1342f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1343f96d37fcSHong Zhang   c->ptap     = ptap;
13440298fd71SBarry Smith   ptap->api   = NULL;
13450298fd71SBarry Smith   ptap->apj   = NULL;
1346f96d37fcSHong Zhang   ptap->merge = merge;
1347d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1348d6ab1dc0SHong Zhang 
1349d6ab1dc0SHong Zhang   *C = Cmpi;
1350d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1351d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
1352d6ab1dc0SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1353d6ab1dc0SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1354d6ab1dc0SHong Zhang   } else {
1355d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1356d6ab1dc0SHong Zhang   }
1357d6ab1dc0SHong Zhang #endif
1358d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1359d6ab1dc0SHong Zhang }
1360d6ab1dc0SHong Zhang 
1361d6ab1dc0SHong Zhang #undef __FUNCT__
1362d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
1363d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C)
1364d6ab1dc0SHong Zhang {
1365d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1366d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1367d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1368d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1369d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1370e745005bSHong Zhang   PetscInt            *adj;
1371e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1372e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1373d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1374*ce94432eSBarry Smith   MPI_Comm            comm;
1375d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1376d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1377d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1378d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1379d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1380d6ab1dc0SHong Zhang   MPI_Status          *status;
1381d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1382d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1383a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1384d6ab1dc0SHong Zhang   Mat                 A_loc;
1385d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1386d6ab1dc0SHong Zhang 
1387d6ab1dc0SHong Zhang   PetscFunctionBegin;
1388*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1389d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1390d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1391d6ab1dc0SHong Zhang 
1392d6ab1dc0SHong Zhang   ptap  = c->ptap;
1393d6ab1dc0SHong Zhang   merge = ptap->merge;
1394d6ab1dc0SHong Zhang 
1395e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1396e745005bSHong Zhang   /*------------------------------------------*/
1397d6ab1dc0SHong Zhang   /* get data from symbolic products */
1398d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1399d6ab1dc0SHong Zhang   ierr   = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1400d6ab1dc0SHong Zhang   ierr   = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1401d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1402d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1403d6ab1dc0SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1404d6ab1dc0SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1405d6ab1dc0SHong Zhang 
1406d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1407d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1408d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1409d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1410d6ab1dc0SHong Zhang   ai    = a_loc->i;
1411d6ab1dc0SHong Zhang   aj    = a_loc->j;
1412d6ab1dc0SHong Zhang 
1413d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1414d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1415d6ab1dc0SHong Zhang     adj = aj + ai[i];
1416d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1417d6ab1dc0SHong Zhang 
1418d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1419e745005bSHong Zhang     /*-------------------------------------------------------------*/
1420d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1421d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1422d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1423d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1424d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1425d6ab1dc0SHong Zhang       row = poJ[j];
1426d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1427d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1428e745005bSHong Zhang       /* perform sparse axpy */
1429e745005bSHong Zhang       nexta  = 0;
1430d6ab1dc0SHong Zhang       valtmp = pA[j];
1431e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1432e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1433e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1434e745005bSHong Zhang           nexta++;
1435d6ab1dc0SHong Zhang         }
1436e745005bSHong Zhang       }
1437e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1438d6ab1dc0SHong Zhang     }
1439d6ab1dc0SHong Zhang 
1440d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1441d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1442d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1443d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1444d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1445d6ab1dc0SHong Zhang       row = pdJ[j];
1446d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1447d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1448e745005bSHong Zhang       /* perform sparse axpy */
1449e745005bSHong Zhang       nexta  = 0;
1450d6ab1dc0SHong Zhang       valtmp = pA[j];
1451e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1452e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1453e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1454e745005bSHong Zhang           nexta++;
1455d6ab1dc0SHong Zhang         }
1456e745005bSHong Zhang       }
1457e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1458d6ab1dc0SHong Zhang     }
1459d6ab1dc0SHong Zhang   }
1460d6ab1dc0SHong Zhang 
1461d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1462d6ab1dc0SHong Zhang   /*------------------------------------*/
1463d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1464d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1465d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1466d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1467d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1468d6ab1dc0SHong Zhang 
1469d6ab1dc0SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1470d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1471d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1472d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1473d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1474d6ab1dc0SHong Zhang     k++;
1475d6ab1dc0SHong Zhang   }
1476d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1477d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1478d6ab1dc0SHong Zhang 
1479d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1480d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1481d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1482d6ab1dc0SHong Zhang 
1483d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1484d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1485d6ab1dc0SHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1486d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1487e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1488d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1489d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1490d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1491d6ab1dc0SHong Zhang   }
1492d6ab1dc0SHong Zhang 
1493d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1494d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1495d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1496d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1497d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1498d6ab1dc0SHong Zhang     /* add received vals into ba */
1499d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1500d6ab1dc0SHong Zhang       /* i-th row */
1501d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1502d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1503d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1504d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1505d6ab1dc0SHong Zhang         nextcj = 0;
1506d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1507d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1508d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1509d6ab1dc0SHong Zhang           }
1510d6ab1dc0SHong Zhang         }
1511d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1512d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1513d6ab1dc0SHong Zhang       }
1514d6ab1dc0SHong Zhang     }
1515d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1516d6ab1dc0SHong Zhang   }
1517d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1518d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1519d6ab1dc0SHong Zhang 
1520d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1521d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1522d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1523d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1524d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1525d6ab1dc0SHong Zhang }
1526d6ab1dc0SHong Zhang 
1527d6ab1dc0SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1528d6ab1dc0SHong Zhang #undef __FUNCT__
1529d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
1530d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C)
1531d6ab1dc0SHong Zhang {
1532d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1533d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1534d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
15350298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1536d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1537d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1538d6ab1dc0SHong Zhang   PetscInt            nnz;
1539d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1540d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1541*ce94432eSBarry Smith   MPI_Comm            comm;
1542d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1543d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1544d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1545d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1546d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1547d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1548d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1549d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1550d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1551d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1552c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1553d6ab1dc0SHong Zhang   PetscScalar         *vals;
1554d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1555d6ab1dc0SHong Zhang 
1556d6ab1dc0SHong Zhang   PetscFunctionBegin;
1557*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1558d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1559d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1560d6ab1dc0SHong 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);
1561d6ab1dc0SHong Zhang   }
1562d6ab1dc0SHong Zhang 
1563d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1564d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1565d6ab1dc0SHong Zhang 
1566d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1567d6ab1dc0SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1568d6ab1dc0SHong Zhang 
1569d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1570d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
15712205254eSKarl Rupp 
1572d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1573d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1574d6ab1dc0SHong Zhang   ai          = a_loc->i;
1575d6ab1dc0SHong Zhang   aj          = a_loc->j;
1576d6ab1dc0SHong Zhang 
1577d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1578d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1579d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1580d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1581d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1582d6ab1dc0SHong Zhang 
1583d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1584d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1585d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1586d6ab1dc0SHong Zhang 
1587d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1588d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1589d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1590d6ab1dc0SHong Zhang   ierr   = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1591d6ab1dc0SHong Zhang   coi[0] = 0;
1592d6ab1dc0SHong Zhang 
1593d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1594d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1595a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1596d6ab1dc0SHong Zhang   current_space = free_space;
159719f0e57aSHong Zhang 
1598d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1599d6ab1dc0SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1600c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1601d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1602d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1603d6ab1dc0SHong Zhang 
1604d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1605d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1606d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1607d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1608d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1609d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1610d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1611d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1612d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1613d6ab1dc0SHong Zhang     }
1614d6ab1dc0SHong Zhang     nnz = lnk[0];
1615d6ab1dc0SHong Zhang 
1616d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1617d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1618d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1619d6ab1dc0SHong Zhang       nspacedouble++;
1620d6ab1dc0SHong Zhang     }
1621d6ab1dc0SHong Zhang 
1622d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1623d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
16242205254eSKarl Rupp 
1625d6ab1dc0SHong Zhang     current_space->array           += nnz;
1626d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1627d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
16282205254eSKarl Rupp 
1629d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1630d6ab1dc0SHong Zhang   }
1631d6ab1dc0SHong Zhang 
1632d6ab1dc0SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1633d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
16342205254eSKarl Rupp 
1635118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1636d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1637d6ab1dc0SHong Zhang 
1638d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1639d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1640d6ab1dc0SHong Zhang   /* determine row ownership */
1641d6ab1dc0SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1642d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
16432205254eSKarl Rupp 
1644d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1645d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
16462205254eSKarl Rupp 
1647d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1648d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1649d6ab1dc0SHong Zhang 
1650d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
1651d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1652d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1653d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
16542205254eSKarl Rupp 
1655d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1656d6ab1dc0SHong Zhang   merge->nsend = 0;
1657d6ab1dc0SHong Zhang 
1658d6ab1dc0SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1659d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1660d6ab1dc0SHong Zhang 
1661d6ab1dc0SHong Zhang   proc = 0;
1662d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1663d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1664d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1665d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1666d6ab1dc0SHong Zhang   }
1667d6ab1dc0SHong Zhang 
1668d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1669d6ab1dc0SHong Zhang   owners_co[0] = 0;
1670d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1671d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1672d6ab1dc0SHong Zhang     if (len_si[proc]) {
1673d6ab1dc0SHong Zhang       merge->nsend++;
1674d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1675d6ab1dc0SHong Zhang       len         += len_si[proc];
1676d6ab1dc0SHong Zhang     }
1677d6ab1dc0SHong Zhang   }
1678d6ab1dc0SHong Zhang 
1679d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
16800298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1681d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1682d6ab1dc0SHong Zhang 
1683d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1684d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1685d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1686d6ab1dc0SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1687d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1688d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1689d6ab1dc0SHong Zhang     i    = owners_co[proc];
1690d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1691d6ab1dc0SHong Zhang     k++;
1692d6ab1dc0SHong Zhang   }
1693d6ab1dc0SHong Zhang 
1694d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1695d6ab1dc0SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1696d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1697c280ed6aSJed Brown     PetscMPIInt icompleted;
1698c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1699d6ab1dc0SHong Zhang   }
1700d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1701d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1702d6ab1dc0SHong Zhang 
1703d6ab1dc0SHong Zhang   /* send and recv coi */
1704d6ab1dc0SHong Zhang   /*-------------------*/
1705d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1706d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1707d6ab1dc0SHong Zhang   ierr   = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1708d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1709d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1710d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1711d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1712d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1713d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1714d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1715d6ab1dc0SHong Zhang     */
1716d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1717d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1718d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1719d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1720d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1721d6ab1dc0SHong Zhang     nrows       = 0;
1722d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1723d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1724d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1725d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1726d6ab1dc0SHong Zhang       nrows++;
1727d6ab1dc0SHong Zhang     }
1728d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1729d6ab1dc0SHong Zhang     k++;
1730d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1731d6ab1dc0SHong Zhang   }
1732d6ab1dc0SHong Zhang   i = merge->nrecv;
1733d6ab1dc0SHong Zhang   while (i--) {
1734c280ed6aSJed Brown     PetscMPIInt icompleted;
1735c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1736d6ab1dc0SHong Zhang   }
1737d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1738d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1739d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1740d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1741d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1742d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1743d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1744d6ab1dc0SHong Zhang 
1745d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1746d6ab1dc0SHong Zhang   /*------------------------------------------*/
1747d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1748d6ab1dc0SHong Zhang   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1749d6ab1dc0SHong Zhang   bi[0] = 0;
1750d6ab1dc0SHong Zhang 
1751d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1752d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1753a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1754d6ab1dc0SHong Zhang   current_space = free_space;
1755d6ab1dc0SHong Zhang 
1756d6ab1dc0SHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1757d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1758d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1759d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1760d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
17612205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1762d6ab1dc0SHong Zhang   }
1763d6ab1dc0SHong Zhang 
1764d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1765d6ab1dc0SHong Zhang   rmax = 0;
1766d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1767d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1768d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1769d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1770d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1771d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1772d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1773d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1774d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1775d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1776d6ab1dc0SHong Zhang     }
1777d6ab1dc0SHong Zhang 
1778d6ab1dc0SHong Zhang     /* add received col data into lnk */
1779d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1780d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1781d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1782d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1783d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1784d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1785d6ab1dc0SHong Zhang       }
1786d6ab1dc0SHong Zhang     }
1787d6ab1dc0SHong Zhang     nnz = lnk[0];
1788d6ab1dc0SHong Zhang 
1789d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1790d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1791d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1792d6ab1dc0SHong Zhang       nspacedouble++;
1793d6ab1dc0SHong Zhang     }
1794d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1795d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1796d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
17972205254eSKarl Rupp 
1798d6ab1dc0SHong Zhang     current_space->array           += nnz;
1799d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1800d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18012205254eSKarl Rupp 
1802d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1803d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1804d6ab1dc0SHong Zhang   }
1805d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1806d6ab1dc0SHong Zhang 
1807d6ab1dc0SHong Zhang   ierr      = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1808d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1809118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1810d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1811d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1812d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1813d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1814d6ab1dc0SHong Zhang 
1815d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1816d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
1817d6ab1dc0SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1818d6ab1dc0SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1819d6ab1dc0SHong Zhang 
1820d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1821d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1822a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,A->cmap->bs);CHKERRQ(ierr);
1823d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1824d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1825d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1826d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1827d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1828d6ab1dc0SHong Zhang     row  = i + rstart;
1829d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1830d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1831d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1832d6ab1dc0SHong Zhang   }
1833d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1834d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1835d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1836d6ab1dc0SHong Zhang 
1837d6ab1dc0SHong Zhang   merge->bi        = bi;
1838d6ab1dc0SHong Zhang   merge->bj        = bj;
1839d6ab1dc0SHong Zhang   merge->coi       = coi;
1840d6ab1dc0SHong Zhang   merge->coj       = coj;
1841d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1842d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1843d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1844d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1845d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1846d6ab1dc0SHong Zhang 
1847d6ab1dc0SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1848d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1849d6ab1dc0SHong Zhang 
1850d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1851d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
18522205254eSKarl Rupp 
1853d6ab1dc0SHong Zhang   c->ptap     = ptap;
18540298fd71SBarry Smith   ptap->api   = NULL;
18550298fd71SBarry Smith   ptap->apj   = NULL;
1856d6ab1dc0SHong Zhang   ptap->merge = merge;
1857d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
18580298fd71SBarry Smith   ptap->apa   = NULL;
1859f96d37fcSHong Zhang 
1860f96d37fcSHong Zhang   *C = Cmpi;
1861f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1862f96d37fcSHong Zhang   if (bi[pn] != 0) {
1863f96d37fcSHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
18641c7d5954SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1865f96d37fcSHong Zhang   } else {
1866f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1867f96d37fcSHong Zhang   }
1868f96d37fcSHong Zhang #endif
1869f96d37fcSHong Zhang   PetscFunctionReturn(0);
1870f96d37fcSHong Zhang }
1871