xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision c5af039c1f267e426bd3017bf1c40568e432060e)
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>
116291f5a6SHong Zhang /*
126291f5a6SHong Zhang #define DEBUG_MATMATMULT
136291f5a6SHong Zhang  */
141c7d5954SHong Zhang /*
151c7d5954SHong Zhang #define DEBUG_MATTrMatMult
161c7d5954SHong Zhang  */
171634b032SHong Zhang 
182499ec78SHong Zhang #undef __FUNCT__
192499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
202499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
212499ec78SHong Zhang {
222499ec78SHong Zhang   PetscErrorCode ierr;
232499ec78SHong Zhang 
242499ec78SHong Zhang   PetscFunctionBegin;
252499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
26598bc09dSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
27a1a4e84aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
28598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
292499ec78SHong Zhang   }
30598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
31598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
32598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
332499ec78SHong Zhang   PetscFunctionReturn(0);
342499ec78SHong Zhang }
352499ec78SHong Zhang 
36f33d1a9aSHong Zhang #undef __FUNCT__
37776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI"
38776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr)
39f33d1a9aSHong Zhang {
40f33d1a9aSHong Zhang   PetscErrorCode       ierr;
419649163dSHong Zhang   Mat_MatMatMultMPI    *mult=(Mat_MatMatMultMPI*)ptr;
42f33d1a9aSHong Zhang 
43f33d1a9aSHong Zhang   PetscFunctionBegin;
446bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr);
456bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr);
466bf464f9SBarry Smith   ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr);
476bf464f9SBarry Smith   ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr);
486bf464f9SBarry Smith   ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr);
496bf464f9SBarry Smith   ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr);
50f33d1a9aSHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
51f33d1a9aSHong Zhang   PetscFunctionReturn(0);
52f33d1a9aSHong Zhang }
53f33d1a9aSHong Zhang 
54598bc09dSHong Zhang #undef __FUNCT__
55a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
56a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
57598bc09dSHong Zhang {
58598bc09dSHong Zhang   PetscErrorCode ierr;
59598bc09dSHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
60598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
61598bc09dSHong Zhang 
62598bc09dSHong Zhang   PetscFunctionBegin;
63598bc09dSHong Zhang   ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
64598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
65a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
66a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
67a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
68a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
69598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
70598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
71598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
72598bc09dSHong Zhang   PetscFunctionReturn(0);
73598bc09dSHong Zhang }
74598bc09dSHong Zhang 
752499ec78SHong Zhang #undef __FUNCT__
76a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32"
77a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A)
782499ec78SHong Zhang {
792499ec78SHong Zhang   PetscErrorCode     ierr;
80776b82aeSLisandro Dalcin   PetscContainer     container;
819649163dSHong Zhang   Mat_MatMatMultMPI  *mult=PETSC_NULL;
822499ec78SHong Zhang 
832499ec78SHong Zhang   PetscFunctionBegin;
849649163dSHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
8529825010SHong Zhang   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist");
86776b82aeSLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
87dce485f0SBarry Smith   A->ops->destroy   = mult->destroy;
88bf0cc555SLisandro Dalcin   A->ops->duplicate = mult->duplicate;
89bf0cc555SLisandro Dalcin   if (A->ops->destroy) {
90992edd73SBarry Smith     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
91bf0cc555SLisandro Dalcin   }
92bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
932499ec78SHong Zhang   PetscFunctionReturn(0);
942499ec78SHong Zhang }
952499ec78SHong Zhang 
962499ec78SHong Zhang #undef __FUNCT__
97a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32"
98a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M)
99b9c92825SBarry Smith {
100abf897f8SHong Zhang   PetscErrorCode     ierr;
101abf897f8SHong Zhang   Mat_MatMatMultMPI  *mult;
102776b82aeSLisandro Dalcin   PetscContainer     container;
103abf897f8SHong Zhang 
104abf897f8SHong Zhang   PetscFunctionBegin;
105abf897f8SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
10629825010SHong Zhang   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist");
107776b82aeSLisandro Dalcin   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
1085c088420SHong Zhang   /* Note: the container is not duplicated, because it requires deep copying of
1095c088420SHong Zhang      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
1105c088420SHong Zhang      These data sets are only used for repeated calling of MatMatMultNumeric().
1115c088420SHong Zhang      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
112dce485f0SBarry Smith   ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr);
113dce485f0SBarry Smith   (*M)->ops->destroy   = mult->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
114dce485f0SBarry Smith   (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */
115abf897f8SHong Zhang   PetscFunctionReturn(0);
116abf897f8SHong Zhang }
117abf897f8SHong Zhang 
118abf897f8SHong Zhang #undef __FUNCT__
119a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
120a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
1214ae0847bSHong Zhang {
1224ae0847bSHong Zhang   PetscErrorCode     ierr;
1234ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
1244ae0847bSHong Zhang   Mat_PtAPMPI        *ptap=a->ptap;
1254ae0847bSHong Zhang 
1264ae0847bSHong Zhang   PetscFunctionBegin;
1274ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
1284ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
1294ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
1304ae0847bSHong Zhang   PetscFunctionReturn(0);
1314ae0847bSHong Zhang }
1324ae0847bSHong Zhang 
1334ae0847bSHong Zhang #undef __FUNCT__
134a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
135a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
136598bc09dSHong Zhang {
137598bc09dSHong Zhang   PetscErrorCode     ierr;
1384ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
139598bc09dSHong Zhang   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
140598bc09dSHong Zhang   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
141598bc09dSHong Zhang   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
142598bc09dSHong Zhang   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
143598bc09dSHong Zhang   Mat_SeqAIJ         *p_loc,*p_oth;
144598bc09dSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
145598bc09dSHong Zhang   PetscScalar        *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
146598bc09dSHong Zhang   PetscInt           cm=C->rmap->n,anz,pnz;
147598bc09dSHong Zhang   Mat_PtAPMPI        *ptap=c->ptap;
14865a57a32SSean Farley   PetscInt           *api,*apj,*apJ,i,j,k,row;
14929825010SHong Zhang   PetscInt           cstart=C->cmap->rstart;
150598bc09dSHong Zhang   PetscInt           cdnz,conz,k0,k1;
1517c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
1527c2e2f26SHong Zhang   PetscMPIInt        rank=a->rank;
1537c2e2f26SHong Zhang #endif
154598bc09dSHong Zhang 
155598bc09dSHong Zhang   PetscFunctionBegin;
1567c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
1577c2e2f26SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank);
1587c2e2f26SHong Zhang #endif
1597c2e2f26SHong Zhang 
160a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
161598bc09dSHong Zhang   /*-----------------------------------------------------*/
162a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
163a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
164a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1656291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
1666291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank);
1676291f5a6SHong Zhang #endif
168598bc09dSHong Zhang 
169598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
170598bc09dSHong Zhang   /*----------------------------------------------------------*/
171598bc09dSHong Zhang   /* get data from symbolic products */
172a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
173a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
174598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
175598bc09dSHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
176598bc09dSHong Zhang 
177598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
178598bc09dSHong Zhang   apa = ptap->apa;
179598bc09dSHong Zhang 
18029825010SHong Zhang   api = ptap->api;
18129825010SHong Zhang   apj = ptap->apj;
182598bc09dSHong Zhang   for (i=0; i<cm; i++) {
183598bc09dSHong Zhang     /* diagonal portion of A */
184598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
185598bc09dSHong Zhang     adj = ad->j + adi[i];
186598bc09dSHong Zhang     ada = ad->a + adi[i];
187598bc09dSHong Zhang     for (j=0; j<anz; j++) {
188598bc09dSHong Zhang       row = adj[j];
189598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
190598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
191598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
192598bc09dSHong Zhang 
193598bc09dSHong Zhang       /* perform dense axpy */
194598bc09dSHong Zhang       valtmp = ada[j];
195598bc09dSHong Zhang       for (k=0; k<pnz; k++){
196598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
197598bc09dSHong Zhang       }
198598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
199598bc09dSHong Zhang     }
200598bc09dSHong Zhang 
201598bc09dSHong Zhang     /* off-diagonal portion of A */
202598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
203598bc09dSHong Zhang     aoj = ao->j + aoi[i];
204598bc09dSHong Zhang     aoa = ao->a + aoi[i];
205598bc09dSHong Zhang     for (j=0; j<anz; j++) {
206598bc09dSHong Zhang       row = aoj[j];
207598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
208598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
209598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
210598bc09dSHong Zhang 
211598bc09dSHong Zhang       /* perform dense axpy */
212598bc09dSHong Zhang       valtmp = aoa[j];
213598bc09dSHong Zhang       for (k=0; k<pnz; k++){
214598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
215598bc09dSHong Zhang       }
216598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
217598bc09dSHong Zhang     }
218598bc09dSHong Zhang 
219598bc09dSHong Zhang     /* set values in C */
220598bc09dSHong Zhang     apJ = apj + api[i];
221598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
222598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
223598bc09dSHong Zhang 
224598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
225598bc09dSHong Zhang     ca = coa + co->i[i];
226598bc09dSHong Zhang     k  = 0;
227598bc09dSHong Zhang     for (k0=0; k0<conz; k0++){
228598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
229598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
230598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
231598bc09dSHong Zhang       k++;
232598bc09dSHong Zhang     }
233598bc09dSHong Zhang 
234598bc09dSHong Zhang     /* diagonal part of C */
235598bc09dSHong Zhang     ca = cda + cd->i[i];
236598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++){
237598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
238598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
239598bc09dSHong Zhang       k++;
240598bc09dSHong Zhang     }
241598bc09dSHong Zhang 
242598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
243598bc09dSHong Zhang     ca = coa + co->i[i];
244598bc09dSHong Zhang     for (; k0<conz; k0++){
245598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
246598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
247598bc09dSHong Zhang       k++;
248598bc09dSHong Zhang     }
249598bc09dSHong Zhang   }
250598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252598bc09dSHong Zhang   PetscFunctionReturn(0);
253598bc09dSHong Zhang }
254598bc09dSHong Zhang 
255598bc09dSHong Zhang #undef __FUNCT__
256a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
257a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
258598bc09dSHong Zhang {
259598bc09dSHong Zhang   PetscErrorCode       ierr;
260598bc09dSHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
261bfcd1a73SHong Zhang   Mat                  Cmpi;
262598bc09dSHong Zhang   Mat_PtAPMPI          *ptap;
263598bc09dSHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
2644ae0847bSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
265bfcd1a73SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2664ae0847bSHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2674ae0847bSHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
268d14fa243SHong Zhang   PetscInt             *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
269bfcd1a73SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
270598bc09dSHong Zhang   PetscBT              lnkbt;
271598bc09dSHong Zhang   PetscScalar          *apa;
272bfcd1a73SHong Zhang   PetscReal            afill;
273cf1a0672SHong Zhang   PetscBool            scalable=PETSC_FALSE;
274f84c536bSHong Zhang   PetscInt             nlnk_max,armax,prmax;
2757c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
2767c2e2f26SHong Zhang   PetscMPIInt          rank=a->rank;
2777c2e2f26SHong Zhang #endif
278598bc09dSHong Zhang 
279598bc09dSHong Zhang   PetscFunctionBegin;
280a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
281cf1a0672SHong 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);
282a1a4e84aSHong Zhang   }
283152983bfSHong Zhang 
284152983bfSHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
285cf1a0672SHong Zhang 
286cf1a0672SHong Zhang     ierr = PetscOptionsBool("-matmatmult_32","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
287cf1a0672SHong Zhang     if (scalable){
288cf1a0672SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(A,P,fill,C);;CHKERRQ(ierr);
289a1a4e84aSHong Zhang       PetscFunctionReturn(0);
290a1a4e84aSHong Zhang     }
291cf1a0672SHong Zhang     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
292cf1a0672SHong Zhang     if (scalable){
29325023028SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr);
2941634b032SHong Zhang       PetscFunctionReturn(0);
2951634b032SHong Zhang     }
296152983bfSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
297a1a4e84aSHong Zhang 
2987c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
2997c2e2f26SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank);
3007c2e2f26SHong Zhang #endif
3017c2e2f26SHong Zhang 
302598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
303598bc09dSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
304598bc09dSHong Zhang 
305598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
306a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
3076291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
3086291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
3096291f5a6SHong Zhang #endif
310598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
311a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
312598bc09dSHong Zhang 
313a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
314a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
315598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
316598bc09dSHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
317598bc09dSHong Zhang 
3186291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
319e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done, Annz %D * P_locnnz %D = %D...\n",rank,ad->i[am]+ao->i[am],p_loc->rmax,(ad->i[am]+ao->i[am])*p_loc->rmax);
3206291f5a6SHong Zhang #endif
3216291f5a6SHong Zhang 
322598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
323598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
324598bc09dSHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
325a1a4e84aSHong Zhang   ptap->api = api;
326598bc09dSHong Zhang   api[0]    = 0;
327598bc09dSHong Zhang 
328598bc09dSHong Zhang   /* create and initialize a linked list */
329f84c536bSHong Zhang   armax = ad->rmax+ao->rmax;
330f84c536bSHong Zhang   prmax = PetscMax(p_loc->rmax,p_oth->rmax);
331f84c536bSHong Zhang   nlnk_max = armax*prmax;
332f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
333f84c536bSHong Zhang #if defined(DEBUG_MATMATMULT)
334f84c536bSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] pN %d; nlnk_max %d; Armax %d+%d=%d; Prmax Max(%d,%d)=%d\n",rank,pN,nlnk_max,ad->rmax,ao->rmax,armax,p_loc->rmax,p_oth->rmax,prmax);
335f84c536bSHong Zhang #endif
3360d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
337598bc09dSHong Zhang 
338bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
339bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
340598bc09dSHong Zhang   current_space = free_space;
341598bc09dSHong Zhang 
342598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
343598bc09dSHong Zhang   for (i=0; i<am; i++) {
344598bc09dSHong Zhang     apnz = 0;
345598bc09dSHong Zhang     /* diagonal portion of A */
346598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
347598bc09dSHong Zhang     for (j=0; j<nzi; j++){
348598bc09dSHong Zhang       row = *adj++;
349598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
350598bc09dSHong Zhang       Jptr  = pj_loc + pi_loc[row];
351598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
3521fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
353598bc09dSHong Zhang     }
354598bc09dSHong Zhang     /* off-diagonal portion of A */
355598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
356598bc09dSHong Zhang     for (j=0; j<nzi; j++){
357598bc09dSHong Zhang       row = *aoj++;
358598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
359598bc09dSHong Zhang       Jptr  = pj_oth + pi_oth[row];
3601fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
361598bc09dSHong Zhang     }
362598bc09dSHong Zhang 
363d14fa243SHong Zhang     apnz     = lnk[0];
364598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
365598bc09dSHong Zhang 
366598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
367598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
368598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
369598bc09dSHong Zhang       nspacedouble++;
370598bc09dSHong Zhang     }
371598bc09dSHong Zhang 
372598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
3731fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
374598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
375598bc09dSHong Zhang     current_space->array           += apnz;
376598bc09dSHong Zhang     current_space->local_used      += apnz;
377598bc09dSHong Zhang     current_space->local_remaining -= apnz;
378598bc09dSHong Zhang   }
379598bc09dSHong Zhang 
380598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
381598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
382a1a4e84aSHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
383a1a4e84aSHong Zhang   apj  = ptap->apj;
384a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
385598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3866291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
3876291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done...\n",rank);
3886291f5a6SHong Zhang #endif
389598bc09dSHong Zhang 
390f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
391f84c536bSHong Zhang   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
392f84c536bSHong Zhang   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
393f84c536bSHong Zhang   ptap->apa = apa;
394f84c536bSHong Zhang #if defined(DEBUG_MATMATMULT)
395f84c536bSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN);
396f84c536bSHong Zhang #endif
397f84c536bSHong Zhang 
398bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
399598bc09dSHong Zhang   /*----------------------------------------------------*/
400bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
401bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
402bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
403bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
404598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
405bfcd1a73SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
406598bc09dSHong Zhang   for (i=0; i<am; i++){
407598bc09dSHong Zhang     row  = i + rstart;
408598bc09dSHong Zhang     apnz = api[i+1] - api[i];
409bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
410598bc09dSHong Zhang     apj += apnz;
411598bc09dSHong Zhang   }
412bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
413bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
414598bc09dSHong Zhang 
415bfcd1a73SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
416bfcd1a73SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
417bfcd1a73SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
418bfcd1a73SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
419bfcd1a73SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
420598bc09dSHong Zhang 
421bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
422bfcd1a73SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
423598bc09dSHong Zhang   c->ptap  = ptap;
424598bc09dSHong Zhang 
425bfcd1a73SHong Zhang   *C = Cmpi;
426bfcd1a73SHong Zhang 
427bfcd1a73SHong Zhang   /* set MatInfo */
428bfcd1a73SHong Zhang   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
429bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
430bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
431bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
432bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
433bfcd1a73SHong Zhang 
434bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
435bfcd1a73SHong Zhang   if (api[am]) {
436bfcd1a73SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
437bfcd1a73SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
438bfcd1a73SHong Zhang   } else {
439bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
440bfcd1a73SHong Zhang   }
441bfcd1a73SHong Zhang #endif
442598bc09dSHong Zhang   PetscFunctionReturn(0);
443598bc09dSHong Zhang }
444598bc09dSHong Zhang 
445a1a4e84aSHong Zhang /* implementation used in PETSc-3.2 */
446a1a4e84aSHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
447598bc09dSHong Zhang #undef __FUNCT__
448cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32"
449cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,Mat C)
450a1a4e84aSHong Zhang {
451a1a4e84aSHong Zhang   PetscErrorCode     ierr;
452a1a4e84aSHong Zhang   Mat                *seq;
453a1a4e84aSHong Zhang   Mat_MatMatMultMPI  *mult;
454a1a4e84aSHong Zhang   PetscContainer     container;
455a1a4e84aSHong Zhang 
456a1a4e84aSHong Zhang   PetscFunctionBegin;
45715f4e0f4SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
45829825010SHong Zhang   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist");
45915f4e0f4SHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
46015f4e0f4SHong Zhang 
46115f4e0f4SHong Zhang   if (mult->skipNumeric){ /* first numeric product is done during symbolic product */
46215f4e0f4SHong Zhang     mult->skipNumeric = PETSC_FALSE;
46315f4e0f4SHong Zhang     PetscFunctionReturn(0);
46415f4e0f4SHong Zhang   }
4657c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
4667c2e2f26SHong Zhang   PetscMPIInt rank;
467fe615159SHong Zhang   MPI_Comm    comm = ((PetscObject)C)->comm;
468fe615159SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
469fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
470cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
4717c2e2f26SHong Zhang #endif
47215f4e0f4SHong Zhang 
473a1a4e84aSHong Zhang   seq = &mult->B_seq;
474a1a4e84aSHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
475a1a4e84aSHong Zhang   mult->B_seq = *seq;
476a1a4e84aSHong Zhang 
477a1a4e84aSHong Zhang   seq = &mult->A_loc;
478a1a4e84aSHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
479a1a4e84aSHong Zhang   mult->A_loc = *seq;
480a1a4e84aSHong Zhang 
48125023028SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
482a1a4e84aSHong Zhang 
483a1a4e84aSHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
484a1a4e84aSHong Zhang   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
485a1a4e84aSHong Zhang   PetscFunctionReturn(0);
486a1a4e84aSHong Zhang }
487a1a4e84aSHong Zhang 
4887c2e2f26SHong Zhang /* numeric product is computed as well */
489a1a4e84aSHong Zhang #undef __FUNCT__
490cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32"
491cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,PetscReal fill,Mat *C)
4922499ec78SHong Zhang {
4932499ec78SHong Zhang   PetscErrorCode     ierr;
4942499ec78SHong Zhang   Mat_MatMatMultMPI  *mult;
495776b82aeSLisandro Dalcin   PetscContainer     container;
496d5821860SHong Zhang   Mat                AB,*seq;
497d5821860SHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
498d5821860SHong Zhang   PetscInt           *idx,i,start,ncols,nzA,nzB,*cmap,imark;
4997c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
500e9e4536cSHong Zhang   MPI_Comm             comm = ((PetscObject)A)->comm;
5017c2e2f26SHong Zhang   PetscMPIInt          rank=a->rank;
5027c2e2f26SHong Zhang #endif
5032499ec78SHong Zhang 
5042499ec78SHong Zhang   PetscFunctionBegin;
5057c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
506cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
5077c2e2f26SHong Zhang #endif
508d0f46423SBarry Smith   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
509cf1a0672SHong Zhang     SETERRQ4(((PetscObject)A)->comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
5102499ec78SHong Zhang   }
511598bc09dSHong Zhang 
5122499ec78SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
5132499ec78SHong Zhang 
514d5821860SHong Zhang   /* get isrowb: nonzero col of A */
515d5821860SHong Zhang   start = A->cmap->rstart;
516d5821860SHong Zhang   cmap  = a->garray;
517d5821860SHong Zhang   nzA   = a->A->cmap->n;
518d5821860SHong Zhang   nzB   = a->B->cmap->n;
519d5821860SHong Zhang   ierr  = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr);
520d5821860SHong Zhang   ncols = 0;
521d5821860SHong Zhang   for (i=0; i<nzB; i++) {  /* row < local row index */
522d5821860SHong Zhang     if (cmap[i] < start) idx[ncols++] = cmap[i];
523d5821860SHong Zhang     else break;
524d5821860SHong Zhang   }
525d5821860SHong Zhang   imark = i;
526d5821860SHong Zhang   for (i=0; i<nzA; i++) idx[ncols++] = start + i;  /* local rows */
527d5821860SHong Zhang   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */
528d5821860SHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr);
529d5821860SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr);
530d5821860SHong Zhang 
531d5821860SHong Zhang   /*  get isrowa: all local rows of A */
532d5821860SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr);
533d5821860SHong Zhang 
534d5821860SHong Zhang   /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */
5352499ec78SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
536d5821860SHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
537d5821860SHong Zhang   mult->B_seq = *seq;
538d5821860SHong Zhang   ierr = PetscFree(seq);CHKERRQ(ierr);
5396291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
5406291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] B_seq is done...\n",rank);
5416291f5a6SHong Zhang #endif
5422499ec78SHong Zhang 
5432499ec78SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
544d5821860SHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
545d5821860SHong Zhang   mult->A_loc = *seq;
546d5821860SHong Zhang   ierr = PetscFree(seq);CHKERRQ(ierr);
5476291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
548e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A_loc is done...\n",rank);
5496291f5a6SHong Zhang #endif
5502499ec78SHong Zhang 
5512499ec78SHong Zhang   /* compute C_seq = A_seq * B_seq */
55225023028SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr);
5536291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
554e9e4536cSHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
555e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Symbolic is done...\n",rank);
556e9e4536cSHong Zhang #endif
55725023028SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
558e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
559e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Numeric is done...\n",rank);
5606291f5a6SHong Zhang #endif
5612499ec78SHong Zhang 
5622499ec78SHong Zhang   /* create mpi matrix C by concatinating C_seq */
5632499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
564*9b8102ccSHong Zhang   ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr);
565*9b8102ccSHong Zhang   ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr);
5666291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
5676291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank);
5686291f5a6SHong Zhang #endif
5692499ec78SHong Zhang 
5702499ec78SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
571776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
572776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
573776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
574d5821860SHong Zhang   ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
575b160f555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
5763f7ae0dbSHong Zhang   mult->skipNumeric  = PETSC_TRUE; /* a numeric product is done here */
577d5821860SHong Zhang   mult->destroy      = AB->ops->destroy;
578d5821860SHong Zhang   mult->duplicate    = AB->ops->duplicate;
579cf1a0672SHong Zhang   AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32;
580a1a4e84aSHong Zhang   AB->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult_32;
581a1a4e84aSHong Zhang   AB->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult_32;
5828cdbd757SHong Zhang   AB->ops->matmult        = MatMatMult_MPIAIJ_MPIAIJ;
583d5821860SHong Zhang   *C                      = AB;
5842499ec78SHong Zhang   PetscFunctionReturn(0);
5852499ec78SHong Zhang }
5862499ec78SHong Zhang 
5879123193aSHong Zhang #undef __FUNCT__
5889123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
5899123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
5909123193aSHong Zhang {
5919123193aSHong Zhang   PetscErrorCode ierr;
5929123193aSHong Zhang 
5939123193aSHong Zhang   PetscFunctionBegin;
5949123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
5959123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
5969123193aSHong Zhang   }
5979123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
5989123193aSHong Zhang   PetscFunctionReturn(0);
5999123193aSHong Zhang }
6009123193aSHong Zhang 
6014324174eSBarry Smith typedef struct {
6024324174eSBarry Smith   Mat         workB;
6034324174eSBarry Smith   PetscScalar *rvalues,*svalues;
6044324174eSBarry Smith   MPI_Request *rwaits,*swaits;
6054324174eSBarry Smith } MPIAIJ_MPIDense;
6064324174eSBarry Smith 
6077af9e4fdSHong Zhang #undef __FUNCT__
6087af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
6094324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
6104324174eSBarry Smith {
6114324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
6124324174eSBarry Smith   PetscErrorCode  ierr;
6134324174eSBarry Smith 
6144324174eSBarry Smith   PetscFunctionBegin;
6156bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
6164324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
6174324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
6184324174eSBarry Smith   PetscFunctionReturn(0);
6194324174eSBarry Smith }
6204324174eSBarry Smith 
6219123193aSHong Zhang #undef __FUNCT__
6229123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
6239123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
6249123193aSHong Zhang {
6259123193aSHong Zhang   PetscErrorCode         ierr;
626f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
627d0f46423SBarry Smith   PetscInt               nz = aij->B->cmap->n;
628bf0cc555SLisandro Dalcin   PetscContainer         container;
6294324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
6304324174eSBarry Smith   VecScatter             ctx = aij->Mvctx;
6314324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
6324324174eSBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
633d0f46423SBarry Smith   PetscInt               m=A->rmap->n,n=B->cmap->n;
6349123193aSHong Zhang 
6359123193aSHong Zhang   PetscFunctionBegin;
636cb2480eaSBarry Smith   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
637d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
638cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
639cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
640cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6418cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense;
642f32d5d43SBarry Smith 
6434324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
644f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
645d0f46423SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
6464324174eSBarry Smith   /* Create work arrays needed */
647d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
648d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
6494324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
6504324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
6514324174eSBarry Smith 
652bf0cc555SLisandro Dalcin   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
653bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
654bf0cc555SLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
655bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
656bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
6579123193aSHong Zhang   PetscFunctionReturn(0);
6589123193aSHong Zhang }
6599123193aSHong Zhang 
6607af9e4fdSHong Zhang #undef __FUNCT__
6617af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
662f32d5d43SBarry Smith /*
6632636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
6642636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
665f32d5d43SBarry Smith */
6664324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
667f32d5d43SBarry Smith {
668f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
669f32d5d43SBarry Smith   PetscErrorCode         ierr;
670f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
671f32d5d43SBarry Smith   VecScatter             ctx = aij->Mvctx;
672f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
673f32d5d43SBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
674f32d5d43SBarry Smith   PetscInt               i,j,k;
675aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
676aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
677f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
6787adad957SLisandro Dalcin   MPI_Comm               comm = ((PetscObject)A)->comm;
679d0f46423SBarry Smith   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
680f32d5d43SBarry Smith   MPI_Status             status;
6814324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
682bf0cc555SLisandro Dalcin   PetscContainer         container;
6834324174eSBarry Smith   Mat                    workB;
684f32d5d43SBarry Smith 
685f32d5d43SBarry Smith   PetscFunctionBegin;
686bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
68729825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
688bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
6894324174eSBarry Smith 
6904324174eSBarry Smith   workB = *outworkB = contents->workB;
691cf1a0672SHong 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);
692f32d5d43SBarry Smith   sindices  = to->indices;
693f32d5d43SBarry Smith   sstarts   = to->starts;
694f32d5d43SBarry Smith   sprocs    = to->procs;
6954324174eSBarry Smith   swaits    = contents->swaits;
6964324174eSBarry Smith   svalues   = contents->svalues;
697f32d5d43SBarry Smith 
698f32d5d43SBarry Smith   rindices  = from->indices;
699f32d5d43SBarry Smith   rstarts   = from->starts;
700f32d5d43SBarry Smith   rprocs    = from->procs;
7014324174eSBarry Smith   rwaits    = contents->rwaits;
7024324174eSBarry Smith   rvalues   = contents->rvalues;
703f32d5d43SBarry Smith 
704f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
705f32d5d43SBarry Smith   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
706f32d5d43SBarry Smith 
707f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
708f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
709f32d5d43SBarry Smith   }
710f32d5d43SBarry Smith 
711f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
712f32d5d43SBarry Smith     /* pack a message at a time */
713f32d5d43SBarry Smith     CHKMEMQ;
714f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
715f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
716f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
7172636ff26SBarry Smith       }
7182636ff26SBarry Smith     }
719f32d5d43SBarry Smith     CHKMEMQ;
720f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
721f32d5d43SBarry Smith   }
7222636ff26SBarry Smith 
723f32d5d43SBarry Smith   nrecvs = from->n;
724f32d5d43SBarry Smith   while (nrecvs) {
725f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
726f32d5d43SBarry Smith     nrecvs--;
727f32d5d43SBarry Smith     /* unpack a message at a time */
728f32d5d43SBarry Smith     CHKMEMQ;
729f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
730f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
731f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
7322636ff26SBarry Smith       }
7332636ff26SBarry Smith     }
734f32d5d43SBarry Smith     CHKMEMQ;
735f32d5d43SBarry Smith   }
736cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
737f32d5d43SBarry Smith 
738f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
739f32d5d43SBarry Smith   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
740f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
741f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
742f32d5d43SBarry Smith   PetscFunctionReturn(0);
743f32d5d43SBarry Smith }
744f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
745f32d5d43SBarry Smith 
7469123193aSHong Zhang #undef __FUNCT__
7479123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
7489123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
7499123193aSHong Zhang {
7509123193aSHong Zhang   PetscErrorCode       ierr;
751f32d5d43SBarry Smith   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
752f32d5d43SBarry Smith   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
753f32d5d43SBarry Smith   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
754f32d5d43SBarry Smith   Mat                  workB;
7559123193aSHong Zhang 
7569123193aSHong Zhang   PetscFunctionBegin;
7579123193aSHong Zhang 
758f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
759f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
760f32d5d43SBarry Smith 
761f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
7624324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
763f32d5d43SBarry Smith 
764f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
765f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
7669123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7679123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7689123193aSHong Zhang   PetscFunctionReturn(0);
7699123193aSHong Zhang }
770cf1a0672SHong Zhang 
7711634b032SHong Zhang #undef __FUNCT__
772cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
773cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
7741634b032SHong Zhang {
7751634b032SHong Zhang   PetscErrorCode     ierr;
776cf1a0672SHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
777cf1a0672SHong Zhang   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
778cf1a0672SHong Zhang   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
779cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
780cf1a0672SHong Zhang   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
781cf1a0672SHong Zhang   Mat_SeqAIJ         *p_loc,*p_oth;
782cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
78329825010SHong Zhang   PetscScalar        *pa_loc,*pa_oth,*pa,valtmp,*ca;
78429825010SHong Zhang   PetscInt           cm=C->rmap->n,anz,pnz;
785cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap=c->ptap;
78629825010SHong Zhang   PetscScalar        *apa_sparse=ptap->apa;
787cf1a0672SHong Zhang   PetscInt           *api,*apj,*apJ,i,j,k,row;
78829825010SHong Zhang   PetscInt           cstart=C->cmap->rstart;
78929825010SHong Zhang   PetscInt           cdnz,conz,k0,k1,nextp;
7901634b032SHong Zhang #if defined(DEBUG_MATMATMULT)
791cf1a0672SHong Zhang   PetscMPIInt        rank=a->rank;
7921634b032SHong Zhang #endif
7931634b032SHong Zhang 
7941634b032SHong Zhang   PetscFunctionBegin;
7951634b032SHong Zhang #if defined(DEBUG_MATMATMULT)
796cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank);
7971634b032SHong Zhang #endif
798cf1a0672SHong Zhang 
799cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
800cf1a0672SHong Zhang   /*-----------------------------------------------------*/
801cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
802cf1a0672SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
803cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
804cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
805cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank);
806cf1a0672SHong Zhang #endif
807cf1a0672SHong Zhang 
808cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
809cf1a0672SHong Zhang   /*----------------------------------------------------------*/
810cf1a0672SHong Zhang   /* get data from symbolic products */
811cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
812cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
813cf1a0672SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
814cf1a0672SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
815cf1a0672SHong Zhang 
81629825010SHong Zhang   api = ptap->api;
81729825010SHong Zhang   apj = ptap->apj;
818cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
81929825010SHong Zhang     apJ = apj + api[i];
82029825010SHong Zhang 
821cf1a0672SHong Zhang     /* diagonal portion of A */
822cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
823cf1a0672SHong Zhang     adj = ad->j + adi[i];
824cf1a0672SHong Zhang     ada = ad->a + adi[i];
825cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
826cf1a0672SHong Zhang       row = adj[j];
827cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
828cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
829cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
83029825010SHong Zhang       /* perform sparse axpy */
831cf1a0672SHong Zhang       valtmp = ada[j];
83229825010SHong Zhang       nextp  = 0;
83329825010SHong Zhang       for (k=0; nextp<pnz; k++) {
83429825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
83529825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
83629825010SHong Zhang         }
8371634b032SHong Zhang       }
838cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
839cf1a0672SHong Zhang     }
8401634b032SHong Zhang 
841cf1a0672SHong Zhang     /* off-diagonal portion of A */
842cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
843cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
844cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
845cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
846cf1a0672SHong Zhang       row = aoj[j];
847cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
848cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
849cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
85029825010SHong Zhang       /* perform sparse axpy */
851cf1a0672SHong Zhang       valtmp = aoa[j];
85229825010SHong Zhang       nextp  = 0;
85329825010SHong Zhang       for (k=0; nextp<pnz; k++) {
85429825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
85529825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
85629825010SHong Zhang         }
857cf1a0672SHong Zhang       }
858cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
859cf1a0672SHong Zhang     }
8601634b032SHong Zhang 
861cf1a0672SHong Zhang     /* set values in C */
862cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
863cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
8641634b032SHong Zhang 
865cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
866cf1a0672SHong Zhang     ca = coa + co->i[i];
867cf1a0672SHong Zhang     k  = 0;
868cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++){
869cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
87029825010SHong Zhang       ca[k0]      = apa_sparse[k];
87129825010SHong Zhang       apa_sparse[k] = 0.0;
872cf1a0672SHong Zhang       k++;
873cf1a0672SHong Zhang     }
8741634b032SHong Zhang 
875cf1a0672SHong Zhang     /* diagonal part of C */
876cf1a0672SHong Zhang     ca = cda + cd->i[i];
877cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++){
87829825010SHong Zhang       ca[k1]      = apa_sparse[k];
87929825010SHong Zhang       apa_sparse[k] = 0.0;
880cf1a0672SHong Zhang       k++;
881cf1a0672SHong Zhang     }
882cf1a0672SHong Zhang 
883cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
884cf1a0672SHong Zhang     ca = coa + co->i[i];
885cf1a0672SHong Zhang     for (; k0<conz; k0++){
88629825010SHong Zhang       ca[k0]      = apa_sparse[k];
88729825010SHong Zhang       apa_sparse[k] = 0.0;
888cf1a0672SHong Zhang       k++;
889cf1a0672SHong Zhang     }
890cf1a0672SHong Zhang   }
891cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
892cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
893cf1a0672SHong Zhang   PetscFunctionReturn(0);
894cf1a0672SHong Zhang }
895cf1a0672SHong Zhang 
896cf1a0672SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */
897cf1a0672SHong Zhang #undef __FUNCT__
898cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
899cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C)
900cf1a0672SHong Zhang {
901cf1a0672SHong Zhang   PetscErrorCode       ierr;
902cf1a0672SHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
903cf1a0672SHong Zhang   Mat                  Cmpi;
904cf1a0672SHong Zhang   Mat_PtAPMPI          *ptap;
905cf1a0672SHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
906cf1a0672SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
907cf1a0672SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
908cf1a0672SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
909cf1a0672SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
910f84c536bSHong Zhang   PetscInt             i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
911cf1a0672SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
912cf1a0672SHong Zhang   PetscInt             nlnk_max,armax,prmax;
913cf1a0672SHong Zhang   PetscReal            afill;
914cf1a0672SHong Zhang   PetscScalar          *apa;
915cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
916cf1a0672SHong Zhang   PetscMPIInt          rank=a->rank;
917cf1a0672SHong Zhang #endif
918cf1a0672SHong Zhang 
919cf1a0672SHong Zhang   PetscFunctionBegin;
920cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
921cf1a0672SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
922cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank);
923cf1a0672SHong Zhang #endif
924cf1a0672SHong Zhang 
925cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
926cf1a0672SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
927cf1a0672SHong Zhang 
928cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
929cf1a0672SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
930cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
931cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
932cf1a0672SHong Zhang #endif
933cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
934cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
935cf1a0672SHong Zhang 
936cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
937cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
938cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
939cf1a0672SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
940cf1a0672SHong Zhang 
941cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
942cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done, Annz %D * P_locnnz %D = %D...\n",rank,ad->i[am]+ao->i[am],p_loc->rmax,(ad->i[am]+ao->i[am])*p_loc->rmax);
943cf1a0672SHong Zhang #endif
944cf1a0672SHong Zhang 
945cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
946cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
947cf1a0672SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
948cf1a0672SHong Zhang   ptap->api = api;
949cf1a0672SHong Zhang   api[0]    = 0;
950cf1a0672SHong Zhang 
951cf1a0672SHong Zhang   /* create and initialize a linked list */
952cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
953cf1a0672SHong Zhang   prmax = PetscMax(p_loc->rmax,p_oth->rmax);
954cf1a0672SHong Zhang   nlnk_max = armax*prmax;
955cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
956cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
957cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] pN %d; nlnk_max %d; Armax %d+%d=%d; Prmax Max(%d,%d)=%d\n",rank,pN,nlnk_max,ad->rmax,ao->rmax,armax,p_loc->rmax,p_oth->rmax,prmax);
958cf1a0672SHong Zhang #endif
959f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
960cf1a0672SHong Zhang 
961cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
962cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
963cf1a0672SHong Zhang   current_space = free_space;
964cf1a0672SHong Zhang 
965cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
966cf1a0672SHong Zhang   for (i=0; i<am; i++) {
967cf1a0672SHong Zhang     apnz = 0;
968cf1a0672SHong Zhang     /* diagonal portion of A */
969cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
970cf1a0672SHong Zhang     for (j=0; j<nzi; j++){
971cf1a0672SHong Zhang       row = *adj++;
972cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
973cf1a0672SHong Zhang       Jptr  = pj_loc + pi_loc[row];
974cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
975f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
976cf1a0672SHong Zhang     }
977cf1a0672SHong Zhang     /* off-diagonal portion of A */
978cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
979cf1a0672SHong Zhang     for (j=0; j<nzi; j++){
980cf1a0672SHong Zhang       row = *aoj++;
981cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
982cf1a0672SHong Zhang       Jptr  = pj_oth + pi_oth[row];
983f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
984cf1a0672SHong Zhang     }
985cf1a0672SHong Zhang 
986f84c536bSHong Zhang     apnz     = *lnk;
987cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
988cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
989cf1a0672SHong Zhang 
990cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
991cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
992cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
993cf1a0672SHong Zhang       nspacedouble++;
994cf1a0672SHong Zhang     }
995cf1a0672SHong Zhang 
996cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
997f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
998cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
999cf1a0672SHong Zhang     current_space->array           += apnz;
1000cf1a0672SHong Zhang     current_space->local_used      += apnz;
1001cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
1002cf1a0672SHong Zhang   }
1003cf1a0672SHong Zhang 
1004cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
1005cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
1006cf1a0672SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
1007cf1a0672SHong Zhang   apj  = ptap->apj;
1008cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
1009f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1010cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
1011cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done..., apnz_max %d\n",rank,apnz_max);
1012cf1a0672SHong Zhang #endif
1013cf1a0672SHong Zhang 
1014cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
1015cf1a0672SHong Zhang   /*----------------------------------------------------*/
1016cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1017cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1018cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1019cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1020cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1021cf1a0672SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1022cf1a0672SHong Zhang 
102329825010SHong Zhang   /* malloc apa for assembly Cmpi */
1024cf1a0672SHong Zhang   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
1025cf1a0672SHong Zhang   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
102629825010SHong Zhang   ptap->apa = apa;
1027cf1a0672SHong Zhang   for (i=0; i<am; i++){
1028cf1a0672SHong Zhang     row  = i + rstart;
1029cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
1030cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
1031cf1a0672SHong Zhang     apj += apnz;
1032cf1a0672SHong Zhang   }
1033cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1034cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1035cf1a0672SHong Zhang 
1036cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
1037cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
1038cf1a0672SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1039cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
1040cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
1041cf1a0672SHong Zhang 
1042cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1043cf1a0672SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
1044cf1a0672SHong Zhang   c->ptap  = ptap;
1045cf1a0672SHong Zhang 
1046cf1a0672SHong Zhang   *C = Cmpi;
1047cf1a0672SHong Zhang 
1048cf1a0672SHong Zhang   /* set MatInfo */
1049cf1a0672SHong Zhang   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
1050cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
1051cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
1052cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
1053cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
1054cf1a0672SHong Zhang 
1055cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
1056cf1a0672SHong Zhang   if (api[am]) {
1057cf1a0672SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1058cf1a0672SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
1059cf1a0672SHong Zhang   } else {
1060cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1061cf1a0672SHong Zhang   }
1062cf1a0672SHong Zhang #endif
10631634b032SHong Zhang   PetscFunctionReturn(0);
10641634b032SHong Zhang }
1065f96d37fcSHong Zhang 
1066f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
1067f96d37fcSHong Zhang #undef __FUNCT__
1068f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
1069c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
1070f96d37fcSHong Zhang {
1071f96d37fcSHong Zhang   PetscErrorCode ierr;
1072f96d37fcSHong Zhang 
1073f96d37fcSHong Zhang   PetscFunctionBegin;
1074f96d37fcSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
1075c5af039cSHong Zhang     ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
1076f96d37fcSHong Zhang   }
1077c5af039cSHong Zhang   ierr = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(P,A,*C);CHKERRQ(ierr);
1078f96d37fcSHong Zhang   PetscFunctionReturn(0);
1079f96d37fcSHong Zhang }
1080f96d37fcSHong Zhang 
1081f96d37fcSHong Zhang #undef __FUNCT__
1082f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
1083c5af039cSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1084f96d37fcSHong Zhang {
1085c5af039cSHong Zhang   PetscErrorCode       ierr;
1086c5af039cSHong Zhang   Mat_Merge_SeqsToMPI  *merge;
1087c5af039cSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1088c5af039cSHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1089c5af039cSHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1090c5af039cSHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
1091c5af039cSHong Zhang   Mat_PtAPMPI          *ptap;
1092c5af039cSHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
1093c5af039cSHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
1094c5af039cSHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
1095c5af039cSHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1096c5af039cSHong Zhang   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1097c5af039cSHong Zhang   MPI_Comm             comm=((PetscObject)C)->comm;
1098c5af039cSHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
1099c5af039cSHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1100c5af039cSHong Zhang   PetscInt             **buf_ri,**buf_rj;
1101c5af039cSHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
1102c5af039cSHong Zhang   MPI_Request          *s_waits,*r_waits;
1103c5af039cSHong Zhang   MPI_Status           *status;
1104c5af039cSHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
1105c5af039cSHong Zhang   PetscInt             *api,*apj,*coi,*coj;
1106c5af039cSHong Zhang   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
1107c5af039cSHong Zhang   PetscInt             sparse_axpy;
1108c5af039cSHong Zhang   Mat                  A_loc;
1109c5af039cSHong Zhang   Mat_SeqAIJ           *a_loc;
1110f96d37fcSHong Zhang 
1111f96d37fcSHong Zhang   PetscFunctionBegin;
1112c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1113c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1114c5af039cSHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1115c5af039cSHong Zhang 
1116c5af039cSHong Zhang   ptap  = c->ptap;
1117c5af039cSHong Zhang   merge = ptap->merge;
1118c5af039cSHong Zhang 
1119c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1120c5af039cSHong Zhang   /*--------------------------------------------------------------*/
1121c5af039cSHong Zhang   /* get data from symbolic products */
1122c5af039cSHong Zhang   coi = merge->coi; coj = merge->coj;
1123c5af039cSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1124c5af039cSHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1125c5af039cSHong Zhang 
1126c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
1127c5af039cSHong Zhang   owners = merge->rowmap->range;
1128c5af039cSHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1129c5af039cSHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1130c5af039cSHong Zhang 
1131c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
1132c5af039cSHong Zhang   A_loc = ptap->A_loc;
1133c5af039cSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1134c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1135c5af039cSHong Zhang   api   = a_loc->i;
1136c5af039cSHong Zhang   apj   = a_loc->j;
1137c5af039cSHong Zhang 
1138c5af039cSHong Zhang   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
1139c5af039cSHong Zhang   ierr = PetscMemzero(apa,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
1140c5af039cSHong Zhang 
1141c5af039cSHong Zhang     for (i=0; i<am; i++) {
1142c5af039cSHong Zhang       /* 2-a) put A[i,:] to dense array apa */
1143c5af039cSHong Zhang       anz = api[i+1] - api[i];
1144c5af039cSHong Zhang       adj = apj + api[i];
1145c5af039cSHong Zhang       ada = a_loc->a + api[i];
1146c5af039cSHong Zhang       for (j=0; j<anz; j++){
1147c5af039cSHong Zhang         apa[adj[j]] = ada[j];
1148c5af039cSHong Zhang       }
1149c5af039cSHong Zhang 
1150c5af039cSHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1151c5af039cSHong Zhang       /*--------------------------------------------------------------*/
1152c5af039cSHong Zhang       /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1153c5af039cSHong Zhang       pnz = po->i[i+1] - po->i[i];
1154c5af039cSHong Zhang       poJ = po->j + po->i[i];
1155c5af039cSHong Zhang       pA  = po->a + po->i[i];
1156c5af039cSHong Zhang       for (j=0; j<pnz; j++){
1157c5af039cSHong Zhang         row = poJ[j];
1158c5af039cSHong Zhang         cnz = coi[row+1] - coi[row];
1159c5af039cSHong Zhang         cj  = coj + coi[row];
1160c5af039cSHong Zhang         ca  = coa + coi[row];
1161c5af039cSHong Zhang         /* perform dense axpy */
1162c5af039cSHong Zhang         valtmp = pA[j];
1163c5af039cSHong Zhang         for (k=0; k<cnz; k++) {
1164c5af039cSHong Zhang           ca[k] += valtmp*apa[cj[k]];
1165c5af039cSHong Zhang         }
1166c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1167c5af039cSHong Zhang       }
1168c5af039cSHong Zhang 
1169c5af039cSHong Zhang       /* put the value into Cd (diagonal part) */
1170c5af039cSHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1171c5af039cSHong Zhang       pdJ = pd->j + pd->i[i];
1172c5af039cSHong Zhang       pA  = pd->a + pd->i[i];
1173c5af039cSHong Zhang       for (j=0; j<pnz; j++){
1174c5af039cSHong Zhang         row = pdJ[j];
1175c5af039cSHong Zhang         cnz = bi[row+1] - bi[row];
1176c5af039cSHong Zhang         cj  = bj + bi[row];
1177c5af039cSHong Zhang         ca  = ba + bi[row];
1178c5af039cSHong Zhang         /* perform dense axpy */
1179c5af039cSHong Zhang         valtmp = pA[j];
1180c5af039cSHong Zhang         for (k=0; k<cnz; k++) {
1181c5af039cSHong Zhang           ca[k] += valtmp*apa[cj[k]];
1182c5af039cSHong Zhang         }
1183c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1184c5af039cSHong Zhang       }
1185c5af039cSHong Zhang 
1186c5af039cSHong Zhang       /* zero the current row of A*P */
1187c5af039cSHong Zhang       apJ = apj + api[i];
1188c5af039cSHong Zhang       for (k=0; k<anz; k++) apa[apJ[k]] = 0.0;
1189c5af039cSHong Zhang     }
1190c5af039cSHong Zhang 
1191c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1192c5af039cSHong Zhang   /*------------------------------------*/
1193c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1194c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1195c5af039cSHong Zhang   len_s  = merge->len_s;
1196c5af039cSHong Zhang   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1197c5af039cSHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1198c5af039cSHong Zhang 
1199c5af039cSHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1200c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++){
1201c5af039cSHong Zhang     if (!len_s[proc]) continue;
1202c5af039cSHong Zhang     i = merge->owners_co[proc];
1203c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1204c5af039cSHong Zhang     k++;
1205c5af039cSHong Zhang   }
1206c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1207c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1208c5af039cSHong Zhang 
1209c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1210c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1211c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1212c5af039cSHong Zhang 
1213c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1214c5af039cSHong Zhang   /*----------------------------------------------------*/
1215c5af039cSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1216c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++){
1217c5af039cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure  -- memory corruption ???*/
1218c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1219c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1220c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
1221c5af039cSHong Zhang   }
1222c5af039cSHong Zhang 
1223c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1224c5af039cSHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
1225c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1226c5af039cSHong Zhang     ba_i = ba + bi[i];
1227c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1228c5af039cSHong Zhang     /* add received vals into ba */
1229c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
1230c5af039cSHong Zhang       /* i-th row */
1231c5af039cSHong Zhang       if (i == *nextrow[k]) {
1232c5af039cSHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
1233c5af039cSHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
1234c5af039cSHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
1235c5af039cSHong Zhang         nextcj = 0;
1236c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++){
1237c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
1238c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1239c5af039cSHong Zhang           }
1240c5af039cSHong Zhang         }
1241c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1242c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1243c5af039cSHong Zhang       }
1244c5af039cSHong Zhang     }
1245c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1246c5af039cSHong Zhang   }
1247c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1248c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1249c5af039cSHong Zhang 
1250c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1251c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1252c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1253c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1254c5af039cSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
1255f96d37fcSHong Zhang   PetscFunctionReturn(0);
1256f96d37fcSHong Zhang }
1257f96d37fcSHong Zhang 
1258c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1259f96d37fcSHong Zhang #undef __FUNCT__
1260f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
1261f96d37fcSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1262f96d37fcSHong Zhang {
1263f96d37fcSHong Zhang   PetscErrorCode       ierr;
1264f96d37fcSHong Zhang   Mat                  Cmpi;
1265f96d37fcSHong Zhang   Mat_PtAPMPI          *ptap;
1266f96d37fcSHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
1267f96d37fcSHong Zhang   Mat_MPIAIJ           *p=(Mat_MPIAIJ*)P->data,*c;
1268f96d37fcSHong Zhang   PetscInt             *pdti,*pdtj,*poti,*potj,*ptJ;
1269f96d37fcSHong Zhang   PetscInt             nnz;
1270f96d37fcSHong Zhang   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1271f96d37fcSHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
1272f96d37fcSHong Zhang   PetscBT              lnkbt;
1273f96d37fcSHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
1274f96d37fcSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1275f96d37fcSHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
1276f96d37fcSHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
1277f96d37fcSHong Zhang   PetscInt             nzi,*bi,*bj;
1278f96d37fcSHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1279f96d37fcSHong Zhang   MPI_Request          *swaits,*rwaits;
1280f96d37fcSHong Zhang   MPI_Status           *sstatus,rstatus;
1281f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI  *merge;
1282f96d37fcSHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
1283f96d37fcSHong Zhang   PetscReal            afill=1.0,afill_tmp;
12841c7d5954SHong Zhang   PetscInt             rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1285f96d37fcSHong Zhang   PetscScalar          *vals;
12861c7d5954SHong Zhang   Mat                  A_loc;
12871c7d5954SHong Zhang   Mat_SeqAIJ           *a_loc;
1288f96d37fcSHong Zhang 
1289f96d37fcSHong Zhang   PetscFunctionBegin;
1290c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1291c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
1292c5af039cSHong 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);
1293c5af039cSHong Zhang   }
1294c5af039cSHong Zhang 
1295f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1296f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12971c7d5954SHong Zhang #if defined(DEBUG_MATTrMatMult)
12981c7d5954SHong Zhang   ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultSymbolic P: %d %d, %d %d; A %d %d, %d %d\n",rank,P->rmap->N,pN,P->rmap->n,P->cmap->n,A->rmap->N,aN,A->rmap->n,A->cmap->n);
12991c7d5954SHong Zhang   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
13001c7d5954SHong Zhang #endif
1301f96d37fcSHong Zhang 
1302f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1303f96d37fcSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1304f96d37fcSHong Zhang 
1305f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1306f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1307c5af039cSHong Zhang   ptap->A_loc = A_loc;
13081c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1309f96d37fcSHong Zhang   api   = a_loc->i;
1310f96d37fcSHong Zhang   apj   = a_loc->j;
1311f96d37fcSHong Zhang 
1312f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1313f96d37fcSHong Zhang   /*----------------------------------------------------*/
1314f96d37fcSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
1315f96d37fcSHong Zhang 
1316f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1317f96d37fcSHong Zhang   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1318f96d37fcSHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1319f96d37fcSHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1320f96d37fcSHong Zhang   coi[0] = 0;
1321f96d37fcSHong Zhang 
1322f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1323f96d37fcSHong Zhang   nnz           = fill*(poti[pon] + api[am]);
1324f96d37fcSHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1325f96d37fcSHong Zhang   current_space = free_space;
13261c7d5954SHong Zhang #if defined(DEBUG_MATTrMatMult)
1327f96d37fcSHong Zhang   ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],api[am]);CHKERRQ(ierr);
1328f96d37fcSHong Zhang   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
13291c7d5954SHong Zhang #endif
1330f96d37fcSHong Zhang   /* create and initialize a linked list */
1331f96d37fcSHong Zhang   nlnk = aN+1;
1332f96d37fcSHong Zhang   ierr = PetscLLCreate(aN,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1333f96d37fcSHong Zhang 
1334f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1335f96d37fcSHong Zhang     nnz = 0;
1336f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1337f96d37fcSHong Zhang     ptJ = potj + poti[i];
1338f96d37fcSHong Zhang     for (j=0; j<pnz; j++){
1339f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1340f96d37fcSHong Zhang       apnz = api[row+1] - api[row];
1341f96d37fcSHong Zhang       Jptr = apj + api[row];
1342f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1343f96d37fcSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1344f96d37fcSHong Zhang       nnz += nlnk;
1345f96d37fcSHong Zhang     }
1346f96d37fcSHong Zhang 
1347f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1348f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1349f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1350f96d37fcSHong Zhang       nspacedouble++;
1351f96d37fcSHong Zhang     }
1352f96d37fcSHong Zhang 
1353f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
1354f96d37fcSHong Zhang     ierr = PetscLLClean(aN,aN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1355f96d37fcSHong Zhang     current_space->array           += nnz;
1356f96d37fcSHong Zhang     current_space->local_used      += nnz;
1357f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
1358f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1359f96d37fcSHong Zhang   }
1360f96d37fcSHong Zhang 
1361f96d37fcSHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1362f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1363f96d37fcSHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]);
1364f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1365f96d37fcSHong Zhang 
1366f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1367f96d37fcSHong Zhang   /*----------------------------------------------*/
1368f96d37fcSHong Zhang   /* determine row ownership */
1369f96d37fcSHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1370f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1371f96d37fcSHong Zhang   merge->rowmap->n = pn;
1372f96d37fcSHong Zhang   merge->rowmap->bs = 1;
1373f96d37fcSHong Zhang   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1374f96d37fcSHong Zhang   owners = merge->rowmap->range;
1375f96d37fcSHong Zhang 
1376f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
1377f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1378f96d37fcSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1379f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1380f96d37fcSHong Zhang   len_s = merge->len_s;
1381f96d37fcSHong Zhang   merge->nsend = 0;
1382f96d37fcSHong Zhang 
1383f96d37fcSHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1384f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1385f96d37fcSHong Zhang 
1386f96d37fcSHong Zhang   proc = 0;
1387f96d37fcSHong Zhang   for (i=0; i<pon; i++){
1388f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1389f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1390f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1391f96d37fcSHong Zhang   }
1392f96d37fcSHong Zhang 
1393f96d37fcSHong Zhang   len   = 0;  /* max length of buf_si[] */
1394f96d37fcSHong Zhang   owners_co[0] = 0;
1395f96d37fcSHong Zhang   for (proc=0; proc<size; proc++){
1396f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1397f96d37fcSHong Zhang     if (len_si[proc]){
1398f96d37fcSHong Zhang       merge->nsend++;
1399f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1400f96d37fcSHong Zhang       len += len_si[proc];
1401f96d37fcSHong Zhang     }
1402f96d37fcSHong Zhang   }
1403f96d37fcSHong Zhang 
1404f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
1405f96d37fcSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1406f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1407f96d37fcSHong Zhang 
1408f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1409f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1410f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1411f96d37fcSHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1412f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++){
1413f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1414f96d37fcSHong Zhang     i = owners_co[proc];
1415f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1416f96d37fcSHong Zhang     k++;
1417f96d37fcSHong Zhang   }
1418f96d37fcSHong Zhang 
1419f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1420f96d37fcSHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1421f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++){
1422f96d37fcSHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
1423f96d37fcSHong Zhang   }
1424f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1425f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1426f96d37fcSHong Zhang 
1427f96d37fcSHong Zhang   /* send and recv coi */
1428f96d37fcSHong Zhang   /*-------------------*/
1429f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1430f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1431f96d37fcSHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1432f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1433f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++){
1434f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1435f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1436f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1437f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1438f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1439f96d37fcSHong Zhang     */
1440f96d37fcSHong Zhang     /*-------------------------------------------*/
1441f96d37fcSHong Zhang     nrows = len_si[proc]/2 - 1;
1442f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1443f96d37fcSHong Zhang     buf_si[0]   = nrows;
1444f96d37fcSHong Zhang     buf_si_i[0] = 0;
1445f96d37fcSHong Zhang     nrows = 0;
1446f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
1447f96d37fcSHong Zhang       nzi = coi[i+1] - coi[i];
1448f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1449f96d37fcSHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
1450f96d37fcSHong Zhang       nrows++;
1451f96d37fcSHong Zhang     }
1452f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1453f96d37fcSHong Zhang     k++;
1454f96d37fcSHong Zhang     buf_si += len_si[proc];
1455f96d37fcSHong Zhang   }
1456f96d37fcSHong Zhang   i = merge->nrecv;
1457f96d37fcSHong Zhang   while (i--) {
1458f96d37fcSHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
1459f96d37fcSHong Zhang   }
1460f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1461f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1462f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1463f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1464f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1465f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1466f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1467f96d37fcSHong Zhang 
1468f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1469f96d37fcSHong Zhang   /*------------------------------------------*/
1470f96d37fcSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1471f96d37fcSHong Zhang 
1472f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1473f96d37fcSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1474f96d37fcSHong Zhang   bi[0] = 0;
1475f96d37fcSHong Zhang 
1476f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1477f96d37fcSHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + api[am]);
1478f96d37fcSHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1479f96d37fcSHong Zhang   current_space = free_space;
1480f96d37fcSHong Zhang 
14811c7d5954SHong Zhang #if defined(DEBUG_MATTrMatMult)
1482f96d37fcSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=%d + %d + %d\n",rank,pdti[pn],poti[pon],api[am]);
14831c7d5954SHong Zhang #endif
1484f96d37fcSHong Zhang 
1485f96d37fcSHong Zhang   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1486f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++){
1487f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1488f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1489f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1490f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
1491f96d37fcSHong Zhang   }
1492f96d37fcSHong Zhang 
14931c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1494f96d37fcSHong Zhang   rmax = 0;
1495f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1496f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1497f96d37fcSHong Zhang     nnz = 0;
1498f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1499f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1500f96d37fcSHong Zhang     for (j=0; j<pnz; j++){
1501f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1502f96d37fcSHong Zhang       apnz = api[row+1] - api[row];
1503f96d37fcSHong Zhang       Jptr = apj + api[row];
1504f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1505f96d37fcSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1506f96d37fcSHong Zhang       nnz += nlnk;
1507f96d37fcSHong Zhang     }
1508f96d37fcSHong Zhang     /* add received col data into lnk */
1509f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
1510f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1511f96d37fcSHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
1512f96d37fcSHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
1513f96d37fcSHong Zhang         ierr = PetscLLAddSorted(nzi,Jptr,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1514f96d37fcSHong Zhang         nnz += nlnk;
1515f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1516f96d37fcSHong Zhang       }
1517f96d37fcSHong Zhang     }
1518f96d37fcSHong Zhang 
1519f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1520f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1521f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1522f96d37fcSHong Zhang       nspacedouble++;
1523f96d37fcSHong Zhang     }
1524f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
1525f96d37fcSHong Zhang     ierr = PetscLLClean(aN,aN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1526f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1527f96d37fcSHong Zhang     current_space->array           += nnz;
1528f96d37fcSHong Zhang     current_space->local_used      += nnz;
1529f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
1530f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1531f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1532f96d37fcSHong Zhang   }
1533f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1534f96d37fcSHong Zhang 
1535f96d37fcSHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1536f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
15371c7d5954SHong Zhang   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + api[am]);
1538f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1539f96d37fcSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
15401c7d5954SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1541f96d37fcSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
1542f96d37fcSHong Zhang 
15431c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
15441c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
15451c7d5954SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
15461c7d5954SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1547f96d37fcSHong Zhang 
1548f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1549f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1550f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1551f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1552f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1553f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1554f96d37fcSHong Zhang   for (i=0; i<pn; i++){
1555f96d37fcSHong Zhang     row = i + rstart;
1556f96d37fcSHong Zhang     nnz = bi[i+1] - bi[i];
1557f96d37fcSHong Zhang     Jptr = bj + bi[i];
1558f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1559f96d37fcSHong Zhang   }
1560f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1561f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1562f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1563f96d37fcSHong Zhang 
1564f96d37fcSHong Zhang   merge->bi            = bi;
1565f96d37fcSHong Zhang   merge->bj            = bj;
1566f96d37fcSHong Zhang   merge->coi           = coi;
1567f96d37fcSHong Zhang   merge->coj           = coj;
1568f96d37fcSHong Zhang   merge->buf_ri        = buf_ri;
1569f96d37fcSHong Zhang   merge->buf_rj        = buf_rj;
1570f96d37fcSHong Zhang   merge->owners_co     = owners_co;
1571f96d37fcSHong Zhang   merge->destroy       = Cmpi->ops->destroy;
1572f96d37fcSHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
1573f96d37fcSHong Zhang 
1574c5af039cSHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
1575f96d37fcSHong Zhang 
1576f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1577f96d37fcSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
1578f96d37fcSHong Zhang   c->ptap        = ptap;
1579c5af039cSHong Zhang   ptap->api      = PETSC_NULL;
1580c5af039cSHong Zhang   ptap->apj      = PETSC_NULL;
1581f96d37fcSHong Zhang   ptap->merge    = merge;
1582f96d37fcSHong Zhang   ptap->apnz_max = ap_rmax;
1583f96d37fcSHong Zhang 
1584f96d37fcSHong Zhang   *C = Cmpi;
1585f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1586f96d37fcSHong Zhang   if (bi[pn] != 0) {
1587f96d37fcSHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
15881c7d5954SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1589f96d37fcSHong Zhang   } else {
1590f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1591f96d37fcSHong Zhang   }
1592f96d37fcSHong Zhang #endif
1593f96d37fcSHong Zhang   PetscFunctionReturn(0);
1594f96d37fcSHong Zhang }
1595