xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision cf1a06720deacc3d28346612b2b7270564520b33)
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  */
141634b032SHong Zhang 
152499ec78SHong Zhang #undef __FUNCT__
162499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
172499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
182499ec78SHong Zhang {
192499ec78SHong Zhang   PetscErrorCode ierr;
202499ec78SHong Zhang 
212499ec78SHong Zhang   PetscFunctionBegin;
222499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
23598bc09dSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
24a1a4e84aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
25598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
262499ec78SHong Zhang   }
27598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
28598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
29598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
302499ec78SHong Zhang   PetscFunctionReturn(0);
312499ec78SHong Zhang }
322499ec78SHong Zhang 
33f33d1a9aSHong Zhang #undef __FUNCT__
34776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI"
35776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr)
36f33d1a9aSHong Zhang {
37f33d1a9aSHong Zhang   PetscErrorCode       ierr;
389649163dSHong Zhang   Mat_MatMatMultMPI    *mult=(Mat_MatMatMultMPI*)ptr;
39f33d1a9aSHong Zhang 
40f33d1a9aSHong Zhang   PetscFunctionBegin;
416bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr);
426bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr);
436bf464f9SBarry Smith   ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr);
446bf464f9SBarry Smith   ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr);
456bf464f9SBarry Smith   ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr);
466bf464f9SBarry Smith   ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr);
47f33d1a9aSHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
48f33d1a9aSHong Zhang   PetscFunctionReturn(0);
49f33d1a9aSHong Zhang }
50f33d1a9aSHong Zhang 
51598bc09dSHong Zhang #undef __FUNCT__
52a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
53a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
54598bc09dSHong Zhang {
55598bc09dSHong Zhang   PetscErrorCode ierr;
56598bc09dSHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
57598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
58598bc09dSHong Zhang 
59598bc09dSHong Zhang   PetscFunctionBegin;
60598bc09dSHong Zhang   ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
61598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
62a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
63a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
64a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
65a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
66598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
67598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
68598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
69598bc09dSHong Zhang   PetscFunctionReturn(0);
70598bc09dSHong Zhang }
71598bc09dSHong Zhang 
722499ec78SHong Zhang #undef __FUNCT__
73a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32"
74a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A)
752499ec78SHong Zhang {
762499ec78SHong Zhang   PetscErrorCode     ierr;
77776b82aeSLisandro Dalcin   PetscContainer     container;
789649163dSHong Zhang   Mat_MatMatMultMPI  *mult=PETSC_NULL;
792499ec78SHong Zhang 
802499ec78SHong Zhang   PetscFunctionBegin;
819649163dSHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
82cf1a0672SHong Zhang   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exit");
83776b82aeSLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
84dce485f0SBarry Smith   A->ops->destroy   = mult->destroy;
85bf0cc555SLisandro Dalcin   A->ops->duplicate = mult->duplicate;
86bf0cc555SLisandro Dalcin   if (A->ops->destroy) {
87992edd73SBarry Smith     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
88bf0cc555SLisandro Dalcin   }
89bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
902499ec78SHong Zhang   PetscFunctionReturn(0);
912499ec78SHong Zhang }
922499ec78SHong Zhang 
932499ec78SHong Zhang #undef __FUNCT__
94a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32"
95a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M)
96b9c92825SBarry Smith {
97abf897f8SHong Zhang   PetscErrorCode     ierr;
98abf897f8SHong Zhang   Mat_MatMatMultMPI  *mult;
99776b82aeSLisandro Dalcin   PetscContainer     container;
100abf897f8SHong Zhang 
101abf897f8SHong Zhang   PetscFunctionBegin;
102abf897f8SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
103cf1a0672SHong Zhang   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exit");
104776b82aeSLisandro Dalcin   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
1055c088420SHong Zhang   /* Note: the container is not duplicated, because it requires deep copying of
1065c088420SHong Zhang      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
1075c088420SHong Zhang      These data sets are only used for repeated calling of MatMatMultNumeric().
1085c088420SHong Zhang      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
109dce485f0SBarry Smith   ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr);
110dce485f0SBarry Smith   (*M)->ops->destroy   = mult->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
111dce485f0SBarry Smith   (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */
112abf897f8SHong Zhang   PetscFunctionReturn(0);
113abf897f8SHong Zhang }
114abf897f8SHong Zhang 
115abf897f8SHong Zhang #undef __FUNCT__
116a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
117a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
1184ae0847bSHong Zhang {
1194ae0847bSHong Zhang   PetscErrorCode     ierr;
1204ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
1214ae0847bSHong Zhang   Mat_PtAPMPI        *ptap=a->ptap;
1224ae0847bSHong Zhang 
1234ae0847bSHong Zhang   PetscFunctionBegin;
1244ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
1254ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
1264ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
1274ae0847bSHong Zhang   PetscFunctionReturn(0);
1284ae0847bSHong Zhang }
1294ae0847bSHong Zhang 
1304ae0847bSHong Zhang #undef __FUNCT__
131a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
132a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
133598bc09dSHong Zhang {
134598bc09dSHong Zhang   PetscErrorCode     ierr;
1354ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
136598bc09dSHong Zhang   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
137598bc09dSHong Zhang   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
138598bc09dSHong Zhang   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
139598bc09dSHong Zhang   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
140598bc09dSHong Zhang   Mat_SeqAIJ         *p_loc,*p_oth;
141598bc09dSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
142598bc09dSHong Zhang   PetscScalar        *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
143598bc09dSHong Zhang   PetscInt           cm=C->rmap->n,anz,pnz;
144598bc09dSHong Zhang   Mat_PtAPMPI        *ptap=c->ptap;
14565a57a32SSean Farley   PetscInt           *api,*apj,*apJ,i,j,k,row;
1464ae0847bSHong Zhang   PetscInt           rstart=C->rmap->rstart,cstart=C->cmap->rstart;
147598bc09dSHong Zhang   PetscInt           cdnz,conz,k0,k1;
1487c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
1497c2e2f26SHong Zhang   PetscMPIInt        rank=a->rank;
1507c2e2f26SHong Zhang #endif
151598bc09dSHong Zhang 
152598bc09dSHong Zhang   PetscFunctionBegin;
1537c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
1547c2e2f26SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank);
1557c2e2f26SHong Zhang #endif
1567c2e2f26SHong Zhang 
157a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
158598bc09dSHong Zhang   /*-----------------------------------------------------*/
159a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
160a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
161a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1626291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
1636291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank);
1646291f5a6SHong Zhang #endif
165598bc09dSHong Zhang 
166598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
167598bc09dSHong Zhang   /*----------------------------------------------------------*/
168598bc09dSHong Zhang   /* get data from symbolic products */
169a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
170a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
171598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
172598bc09dSHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
173598bc09dSHong Zhang 
174598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
175598bc09dSHong Zhang   apa = ptap->apa;
176598bc09dSHong Zhang 
177598bc09dSHong Zhang   for (i=0; i<cm; i++) {
178598bc09dSHong Zhang     /* diagonal portion of A */
179598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
180598bc09dSHong Zhang     adj = ad->j + adi[i];
181598bc09dSHong Zhang     ada = ad->a + adi[i];
182598bc09dSHong Zhang     for (j=0; j<anz; j++) {
183598bc09dSHong Zhang       row = adj[j];
184598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
185598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
186598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
187598bc09dSHong Zhang 
188598bc09dSHong Zhang       /* perform dense axpy */
189598bc09dSHong Zhang       valtmp = ada[j];
190598bc09dSHong Zhang       for (k=0; k<pnz; k++){
191598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
192598bc09dSHong Zhang       }
193598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
194598bc09dSHong Zhang     }
195598bc09dSHong Zhang 
196598bc09dSHong Zhang     /* off-diagonal portion of A */
197598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
198598bc09dSHong Zhang     aoj = ao->j + aoi[i];
199598bc09dSHong Zhang     aoa = ao->a + aoi[i];
200598bc09dSHong Zhang     for (j=0; j<anz; j++) {
201598bc09dSHong Zhang       row = aoj[j];
202598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
203598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
204598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
205598bc09dSHong Zhang 
206598bc09dSHong Zhang       /* perform dense axpy */
207598bc09dSHong Zhang       valtmp = aoa[j];
208598bc09dSHong Zhang       for (k=0; k<pnz; k++){
209598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
210598bc09dSHong Zhang       }
211598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
212598bc09dSHong Zhang     }
213598bc09dSHong Zhang 
214598bc09dSHong Zhang     /* set values in C */
215598bc09dSHong Zhang     row = rstart + i;
216a1a4e84aSHong Zhang     api = ptap->api;
217a1a4e84aSHong Zhang     apj = ptap->apj;
218598bc09dSHong Zhang     apJ = apj + api[i];
219598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
220598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
221598bc09dSHong Zhang 
222598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
223598bc09dSHong Zhang     ca = coa + co->i[i];
224598bc09dSHong Zhang     k  = 0;
225598bc09dSHong Zhang     for (k0=0; k0<conz; k0++){
226598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
227598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
228598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
229598bc09dSHong Zhang       k++;
230598bc09dSHong Zhang     }
231598bc09dSHong Zhang 
232598bc09dSHong Zhang     /* diagonal part of C */
233598bc09dSHong Zhang     ca = cda + cd->i[i];
234598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++){
235598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
236598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
237598bc09dSHong Zhang       k++;
238598bc09dSHong Zhang     }
239598bc09dSHong Zhang 
240598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
241598bc09dSHong Zhang     ca = coa + co->i[i];
242598bc09dSHong Zhang     for (; k0<conz; k0++){
243598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
244598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
245598bc09dSHong Zhang       k++;
246598bc09dSHong Zhang     }
247598bc09dSHong Zhang   }
248598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2506291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
2516291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank);
2526291f5a6SHong Zhang #endif
253598bc09dSHong Zhang   PetscFunctionReturn(0);
254598bc09dSHong Zhang }
255598bc09dSHong Zhang 
256598bc09dSHong Zhang #undef __FUNCT__
257a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
258a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
259598bc09dSHong Zhang {
260598bc09dSHong Zhang   PetscErrorCode       ierr;
261598bc09dSHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
262bfcd1a73SHong Zhang   Mat                  Cmpi;
263598bc09dSHong Zhang   Mat_PtAPMPI          *ptap;
264598bc09dSHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
2654ae0847bSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
266bfcd1a73SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2674ae0847bSHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2684ae0847bSHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
269bfcd1a73SHong Zhang   PetscInt             nlnk,*lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
270bfcd1a73SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
271598bc09dSHong Zhang   PetscBT              lnkbt;
272598bc09dSHong Zhang   PetscScalar          *apa;
273bfcd1a73SHong Zhang   PetscReal            afill;
274cf1a0672SHong Zhang   PetscBool            scalable=PETSC_FALSE;
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   /* malloc apa to store dense row A[i,:]*P */
306bfcd1a73SHong Zhang   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
307bfcd1a73SHong Zhang   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
308598bc09dSHong Zhang   ptap->apa = apa;
3096291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
310e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN);
3116291f5a6SHong Zhang #endif
312598bc09dSHong Zhang 
313598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
314a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
3156291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
3166291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
3176291f5a6SHong Zhang #endif
318598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
319a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
320598bc09dSHong Zhang 
321a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
322a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
323598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
324598bc09dSHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
325598bc09dSHong Zhang 
3266291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
327e9e4536cSHong 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);
3286291f5a6SHong Zhang #endif
3296291f5a6SHong Zhang 
330598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
331598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
332598bc09dSHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
333a1a4e84aSHong Zhang   ptap->api = api;
334598bc09dSHong Zhang   api[0]    = 0;
335598bc09dSHong Zhang 
336598bc09dSHong Zhang   /* create and initialize a linked list */
337598bc09dSHong Zhang   nlnk = pN+1;
338598bc09dSHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
339598bc09dSHong Zhang 
340bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
341bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
342598bc09dSHong Zhang   current_space = free_space;
343598bc09dSHong Zhang 
344598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
345598bc09dSHong Zhang   for (i=0; i<am; i++) {
346598bc09dSHong Zhang     apnz = 0;
347598bc09dSHong Zhang     /* diagonal portion of A */
348598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
349598bc09dSHong Zhang     for (j=0; j<nzi; j++){
350598bc09dSHong Zhang       row = *adj++;
351598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
352598bc09dSHong Zhang       Jptr  = pj_loc + pi_loc[row];
353598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
354dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
355598bc09dSHong Zhang       apnz += nlnk;
356598bc09dSHong Zhang     }
357598bc09dSHong Zhang     /* off-diagonal portion of A */
358598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
359598bc09dSHong Zhang     for (j=0; j<nzi; j++){
360598bc09dSHong Zhang       row = *aoj++;
361598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
362598bc09dSHong Zhang       Jptr  = pj_oth + pi_oth[row];
363dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
364598bc09dSHong Zhang       apnz += nlnk;
365598bc09dSHong Zhang     }
366598bc09dSHong Zhang 
367598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
368598bc09dSHong Zhang 
369598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
370598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
371598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
372598bc09dSHong Zhang       nspacedouble++;
373598bc09dSHong Zhang     }
374598bc09dSHong Zhang 
375598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
376598bc09dSHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
377598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
378598bc09dSHong Zhang     current_space->array           += apnz;
379598bc09dSHong Zhang     current_space->local_used      += apnz;
380598bc09dSHong Zhang     current_space->local_remaining -= apnz;
381598bc09dSHong Zhang   }
382598bc09dSHong Zhang 
383598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
384598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
385a1a4e84aSHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
386a1a4e84aSHong Zhang   apj  = ptap->apj;
387a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
388598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3896291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
3906291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done...\n",rank);
3916291f5a6SHong Zhang #endif
392598bc09dSHong Zhang 
393bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
394598bc09dSHong Zhang   /*----------------------------------------------------*/
395bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
396bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
397bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
398bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
399598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
400bfcd1a73SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
401598bc09dSHong Zhang   for (i=0; i<am; i++){
402598bc09dSHong Zhang     row  = i + rstart;
403598bc09dSHong Zhang     apnz = api[i+1] - api[i];
404bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
405598bc09dSHong Zhang     apj += apnz;
406598bc09dSHong Zhang   }
407bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
408bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409598bc09dSHong Zhang 
410bfcd1a73SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
411bfcd1a73SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
412bfcd1a73SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
413bfcd1a73SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
414bfcd1a73SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
415598bc09dSHong Zhang 
416bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
417bfcd1a73SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
418598bc09dSHong Zhang   c->ptap  = ptap;
419598bc09dSHong Zhang 
420bfcd1a73SHong Zhang   *C = Cmpi;
421bfcd1a73SHong Zhang 
422bfcd1a73SHong Zhang   /* set MatInfo */
423bfcd1a73SHong Zhang   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
424bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
425bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
426bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
427bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
428bfcd1a73SHong Zhang 
429bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
430bfcd1a73SHong Zhang   if (api[am]) {
431bfcd1a73SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
432bfcd1a73SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
433bfcd1a73SHong Zhang   } else {
434bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
435bfcd1a73SHong Zhang   }
436bfcd1a73SHong Zhang #endif
4376291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
4386291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank);
4396291f5a6SHong Zhang #endif
440598bc09dSHong Zhang   PetscFunctionReturn(0);
441598bc09dSHong Zhang }
442598bc09dSHong Zhang 
443a1a4e84aSHong Zhang /* implementation used in PETSc-3.2 */
444a1a4e84aSHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
445598bc09dSHong Zhang #undef __FUNCT__
446cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32"
447cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,Mat C)
448a1a4e84aSHong Zhang {
449a1a4e84aSHong Zhang   PetscErrorCode     ierr;
450a1a4e84aSHong Zhang   Mat                *seq;
451a1a4e84aSHong Zhang   Mat_MatMatMultMPI  *mult;
452a1a4e84aSHong Zhang   PetscContainer     container;
453a1a4e84aSHong Zhang 
454a1a4e84aSHong Zhang   PetscFunctionBegin;
45515f4e0f4SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
456cf1a0672SHong Zhang   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exit");
45715f4e0f4SHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
45815f4e0f4SHong Zhang 
45915f4e0f4SHong Zhang   if (mult->skipNumeric){ /* first numeric product is done during symbolic product */
46015f4e0f4SHong Zhang     mult->skipNumeric = PETSC_FALSE;
46115f4e0f4SHong Zhang     PetscFunctionReturn(0);
46215f4e0f4SHong Zhang   }
4637c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
4647c2e2f26SHong Zhang   PetscMPIInt rank;
465fe615159SHong Zhang   MPI_Comm    comm = ((PetscObject)C)->comm;
466fe615159SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
467fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
468cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
4697c2e2f26SHong Zhang #endif
47015f4e0f4SHong Zhang 
471a1a4e84aSHong Zhang   seq = &mult->B_seq;
472a1a4e84aSHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
473a1a4e84aSHong Zhang   mult->B_seq = *seq;
474a1a4e84aSHong Zhang 
475a1a4e84aSHong Zhang   seq = &mult->A_loc;
476a1a4e84aSHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
477a1a4e84aSHong Zhang   mult->A_loc = *seq;
478a1a4e84aSHong Zhang 
47925023028SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
480a1a4e84aSHong Zhang 
481a1a4e84aSHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
482a1a4e84aSHong Zhang   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
483fe615159SHong Zhang #if defined(DEBUG_MATMATMULT)
484fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
485cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
486fe615159SHong Zhang #endif
487a1a4e84aSHong Zhang   PetscFunctionReturn(0);
488a1a4e84aSHong Zhang }
489a1a4e84aSHong Zhang 
4907c2e2f26SHong Zhang /* numeric product is computed as well */
491a1a4e84aSHong Zhang #undef __FUNCT__
492cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32"
493cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,PetscReal fill,Mat *C)
4942499ec78SHong Zhang {
4952499ec78SHong Zhang   PetscErrorCode     ierr;
4962499ec78SHong Zhang   Mat_MatMatMultMPI  *mult;
497776b82aeSLisandro Dalcin   PetscContainer     container;
498d5821860SHong Zhang   Mat                AB,*seq;
499d5821860SHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
500d5821860SHong Zhang   PetscInt           *idx,i,start,ncols,nzA,nzB,*cmap,imark;
5017c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
502e9e4536cSHong Zhang   MPI_Comm             comm = ((PetscObject)A)->comm;
5037c2e2f26SHong Zhang   PetscMPIInt          rank=a->rank;
5047c2e2f26SHong Zhang #endif
5052499ec78SHong Zhang 
5062499ec78SHong Zhang   PetscFunctionBegin;
5077c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT)
508cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
5097c2e2f26SHong Zhang #endif
510d0f46423SBarry Smith   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
511cf1a0672SHong 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);
5122499ec78SHong Zhang   }
513598bc09dSHong Zhang 
5142499ec78SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
5152499ec78SHong Zhang 
516d5821860SHong Zhang   /* get isrowb: nonzero col of A */
517d5821860SHong Zhang   start = A->cmap->rstart;
518d5821860SHong Zhang   cmap  = a->garray;
519d5821860SHong Zhang   nzA   = a->A->cmap->n;
520d5821860SHong Zhang   nzB   = a->B->cmap->n;
521d5821860SHong Zhang   ierr  = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr);
522d5821860SHong Zhang   ncols = 0;
523d5821860SHong Zhang   for (i=0; i<nzB; i++) {  /* row < local row index */
524d5821860SHong Zhang     if (cmap[i] < start) idx[ncols++] = cmap[i];
525d5821860SHong Zhang     else break;
526d5821860SHong Zhang   }
527d5821860SHong Zhang   imark = i;
528d5821860SHong Zhang   for (i=0; i<nzA; i++) idx[ncols++] = start + i;  /* local rows */
529d5821860SHong Zhang   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */
530d5821860SHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr);
531d5821860SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr);
532d5821860SHong Zhang 
533d5821860SHong Zhang   /*  get isrowa: all local rows of A */
534d5821860SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr);
535d5821860SHong Zhang 
536d5821860SHong Zhang   /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */
5372499ec78SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
538d5821860SHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
539d5821860SHong Zhang   mult->B_seq = *seq;
540d5821860SHong Zhang   ierr = PetscFree(seq);CHKERRQ(ierr);
5416291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
5426291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] B_seq is done...\n",rank);
5436291f5a6SHong Zhang #endif
5442499ec78SHong Zhang 
5452499ec78SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
546d5821860SHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
547d5821860SHong Zhang   mult->A_loc = *seq;
548d5821860SHong Zhang   ierr = PetscFree(seq);CHKERRQ(ierr);
5496291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
550e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A_loc is done...\n",rank);
5516291f5a6SHong Zhang #endif
5522499ec78SHong Zhang 
5532499ec78SHong Zhang   /* compute C_seq = A_seq * B_seq */
55425023028SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr);
5556291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
556e9e4536cSHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
557e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Symbolic is done...\n",rank);
558e9e4536cSHong Zhang #endif
55925023028SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
560e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT)
561e9e4536cSHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Numeric is done...\n",rank);
5626291f5a6SHong Zhang #endif
5632499ec78SHong Zhang 
5642499ec78SHong Zhang   /* create mpi matrix C by concatinating C_seq */
5652499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
566*9b8102ccSHong Zhang   ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr);
567*9b8102ccSHong Zhang   ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr);
5686291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
5696291f5a6SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank);
5706291f5a6SHong Zhang #endif
5712499ec78SHong Zhang 
5722499ec78SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
573776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
574776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
575776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
576d5821860SHong Zhang   ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
577b160f555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
5783f7ae0dbSHong Zhang   mult->skipNumeric  = PETSC_TRUE; /* a numeric product is done here */
579d5821860SHong Zhang   mult->destroy      = AB->ops->destroy;
580d5821860SHong Zhang   mult->duplicate    = AB->ops->duplicate;
581cf1a0672SHong Zhang   AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32;
582a1a4e84aSHong Zhang   AB->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult_32;
583a1a4e84aSHong Zhang   AB->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult_32;
5848cdbd757SHong Zhang   AB->ops->matmult        = MatMatMult_MPIAIJ_MPIAIJ;
585d5821860SHong Zhang   *C                      = AB;
5866291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT)
587cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
5886291f5a6SHong Zhang #endif
5892499ec78SHong Zhang   PetscFunctionReturn(0);
5902499ec78SHong Zhang }
5912499ec78SHong Zhang 
5929123193aSHong Zhang #undef __FUNCT__
5939123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
5949123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
5959123193aSHong Zhang {
5969123193aSHong Zhang   PetscErrorCode ierr;
5979123193aSHong Zhang 
5989123193aSHong Zhang   PetscFunctionBegin;
5999123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6009123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
6019123193aSHong Zhang   }
6029123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
6039123193aSHong Zhang   PetscFunctionReturn(0);
6049123193aSHong Zhang }
6059123193aSHong Zhang 
6064324174eSBarry Smith typedef struct {
6074324174eSBarry Smith   Mat         workB;
6084324174eSBarry Smith   PetscScalar *rvalues,*svalues;
6094324174eSBarry Smith   MPI_Request *rwaits,*swaits;
6104324174eSBarry Smith } MPIAIJ_MPIDense;
6114324174eSBarry Smith 
6127af9e4fdSHong Zhang #undef __FUNCT__
6137af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
6144324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
6154324174eSBarry Smith {
6164324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
6174324174eSBarry Smith   PetscErrorCode  ierr;
6184324174eSBarry Smith 
6194324174eSBarry Smith   PetscFunctionBegin;
6206bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
6214324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
6224324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
6234324174eSBarry Smith   PetscFunctionReturn(0);
6244324174eSBarry Smith }
6254324174eSBarry Smith 
6269123193aSHong Zhang #undef __FUNCT__
6279123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
6289123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
6299123193aSHong Zhang {
6309123193aSHong Zhang   PetscErrorCode         ierr;
631f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
632d0f46423SBarry Smith   PetscInt               nz = aij->B->cmap->n;
633bf0cc555SLisandro Dalcin   PetscContainer         container;
6344324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
6354324174eSBarry Smith   VecScatter             ctx = aij->Mvctx;
6364324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
6374324174eSBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
638d0f46423SBarry Smith   PetscInt               m=A->rmap->n,n=B->cmap->n;
6399123193aSHong Zhang 
6409123193aSHong Zhang   PetscFunctionBegin;
641cb2480eaSBarry Smith   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
642d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
643cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
644cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
645cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6468cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense;
647f32d5d43SBarry Smith 
6484324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
649f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
650d0f46423SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
6514324174eSBarry Smith   /* Create work arrays needed */
652d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
653d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
6544324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
6554324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
6564324174eSBarry Smith 
657bf0cc555SLisandro Dalcin   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
658bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
659bf0cc555SLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
660bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
661bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
6629123193aSHong Zhang   PetscFunctionReturn(0);
6639123193aSHong Zhang }
6649123193aSHong Zhang 
6657af9e4fdSHong Zhang #undef __FUNCT__
6667af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
667f32d5d43SBarry Smith /*
6682636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
6692636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
670f32d5d43SBarry Smith */
6714324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
672f32d5d43SBarry Smith {
673f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
674f32d5d43SBarry Smith   PetscErrorCode         ierr;
675f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
676f32d5d43SBarry Smith   VecScatter             ctx = aij->Mvctx;
677f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
678f32d5d43SBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
679f32d5d43SBarry Smith   PetscInt               i,j,k;
680aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
681aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
682f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
6837adad957SLisandro Dalcin   MPI_Comm               comm = ((PetscObject)A)->comm;
684d0f46423SBarry Smith   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
685f32d5d43SBarry Smith   MPI_Status             status;
6864324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
687bf0cc555SLisandro Dalcin   PetscContainer         container;
6884324174eSBarry Smith   Mat                    workB;
689f32d5d43SBarry Smith 
690f32d5d43SBarry Smith   PetscFunctionBegin;
691bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
692cf1a0672SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exit");
693bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
6944324174eSBarry Smith 
6954324174eSBarry Smith   workB = *outworkB = contents->workB;
696cf1a0672SHong 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);
697f32d5d43SBarry Smith   sindices  = to->indices;
698f32d5d43SBarry Smith   sstarts   = to->starts;
699f32d5d43SBarry Smith   sprocs    = to->procs;
7004324174eSBarry Smith   swaits    = contents->swaits;
7014324174eSBarry Smith   svalues   = contents->svalues;
702f32d5d43SBarry Smith 
703f32d5d43SBarry Smith   rindices  = from->indices;
704f32d5d43SBarry Smith   rstarts   = from->starts;
705f32d5d43SBarry Smith   rprocs    = from->procs;
7064324174eSBarry Smith   rwaits    = contents->rwaits;
7074324174eSBarry Smith   rvalues   = contents->rvalues;
708f32d5d43SBarry Smith 
709f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
710f32d5d43SBarry Smith   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
711f32d5d43SBarry Smith 
712f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
713f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
714f32d5d43SBarry Smith   }
715f32d5d43SBarry Smith 
716f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
717f32d5d43SBarry Smith     /* pack a message at a time */
718f32d5d43SBarry Smith     CHKMEMQ;
719f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
720f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
721f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
7222636ff26SBarry Smith       }
7232636ff26SBarry Smith     }
724f32d5d43SBarry Smith     CHKMEMQ;
725f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
726f32d5d43SBarry Smith   }
7272636ff26SBarry Smith 
728f32d5d43SBarry Smith   nrecvs = from->n;
729f32d5d43SBarry Smith   while (nrecvs) {
730f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
731f32d5d43SBarry Smith     nrecvs--;
732f32d5d43SBarry Smith     /* unpack a message at a time */
733f32d5d43SBarry Smith     CHKMEMQ;
734f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
735f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
736f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
7372636ff26SBarry Smith       }
7382636ff26SBarry Smith     }
739f32d5d43SBarry Smith     CHKMEMQ;
740f32d5d43SBarry Smith   }
741cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
742f32d5d43SBarry Smith 
743f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
744f32d5d43SBarry Smith   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
745f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
746f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
747f32d5d43SBarry Smith   PetscFunctionReturn(0);
748f32d5d43SBarry Smith }
749f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
750f32d5d43SBarry Smith 
7519123193aSHong Zhang #undef __FUNCT__
7529123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
7539123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
7549123193aSHong Zhang {
7559123193aSHong Zhang   PetscErrorCode       ierr;
756f32d5d43SBarry Smith   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
757f32d5d43SBarry Smith   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
758f32d5d43SBarry Smith   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
759f32d5d43SBarry Smith   Mat                  workB;
7609123193aSHong Zhang 
7619123193aSHong Zhang   PetscFunctionBegin;
7629123193aSHong Zhang 
763f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
764f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
765f32d5d43SBarry Smith 
766f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
7674324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
768f32d5d43SBarry Smith 
769f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
770f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
7719123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7729123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7739123193aSHong Zhang   PetscFunctionReturn(0);
7749123193aSHong Zhang }
775cf1a0672SHong Zhang 
7761634b032SHong Zhang #undef __FUNCT__
777cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
778cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
7791634b032SHong Zhang {
7801634b032SHong Zhang   PetscErrorCode     ierr;
781cf1a0672SHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
782cf1a0672SHong Zhang   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
783cf1a0672SHong Zhang   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
784cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
785cf1a0672SHong Zhang   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
786cf1a0672SHong Zhang   Mat_SeqAIJ         *p_loc,*p_oth;
787cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
788cf1a0672SHong Zhang   PetscScalar        *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
789cf1a0672SHong Zhang   PetscInt           cm=C->rmap->n,anz,pnz,pN=P->cmap->N;
790cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap=c->ptap;
791cf1a0672SHong Zhang   PetscInt           *api,*apj,*apJ,i,j,k,row;
792cf1a0672SHong Zhang   PetscInt           rstart=C->rmap->rstart,cstart=C->cmap->rstart;
793cf1a0672SHong Zhang   PetscInt           cdnz,conz,k0,k1;
7941634b032SHong Zhang #if defined(DEBUG_MATMATMULT)
795cf1a0672SHong Zhang   PetscMPIInt        rank=a->rank;
7961634b032SHong Zhang #endif
7971634b032SHong Zhang 
7981634b032SHong Zhang   PetscFunctionBegin;
7991634b032SHong Zhang #if defined(DEBUG_MATMATMULT)
800cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank);
8011634b032SHong Zhang #endif
802cf1a0672SHong Zhang 
803cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
804cf1a0672SHong Zhang   /*-----------------------------------------------------*/
805cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
806cf1a0672SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
807cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
808cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
809cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank);
810cf1a0672SHong Zhang #endif
811cf1a0672SHong Zhang 
812cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
813cf1a0672SHong Zhang   /*----------------------------------------------------------*/
814cf1a0672SHong Zhang   /* get data from symbolic products */
815cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
816cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
817cf1a0672SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
818cf1a0672SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
819cf1a0672SHong Zhang 
820cf1a0672SHong Zhang   /* get apa for storing dense row A[i,:]*P */
821cf1a0672SHong Zhang   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
822cf1a0672SHong Zhang   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
823cf1a0672SHong Zhang 
824cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
825cf1a0672SHong Zhang     /* diagonal portion of A */
826cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
827cf1a0672SHong Zhang     adj = ad->j + adi[i];
828cf1a0672SHong Zhang     ada = ad->a + adi[i];
829cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
830cf1a0672SHong Zhang       row = adj[j];
831cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
832cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
833cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
834cf1a0672SHong Zhang 
835cf1a0672SHong Zhang       /* perform dense axpy */
836cf1a0672SHong Zhang       valtmp = ada[j];
837cf1a0672SHong Zhang       for (k=0; k<pnz; k++){
838cf1a0672SHong Zhang         apa[pj[k]] += valtmp*pa[k];
8391634b032SHong Zhang       }
840cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
841cf1a0672SHong Zhang     }
8421634b032SHong Zhang 
843cf1a0672SHong Zhang     /* off-diagonal portion of A */
844cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
845cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
846cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
847cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
848cf1a0672SHong Zhang       row = aoj[j];
849cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
850cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
851cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
8521634b032SHong Zhang 
853cf1a0672SHong Zhang       /* perform dense axpy */
854cf1a0672SHong Zhang       valtmp = aoa[j];
855cf1a0672SHong Zhang       for (k=0; k<pnz; k++){
856cf1a0672SHong Zhang         apa[pj[k]] += valtmp*pa[k];
857cf1a0672SHong Zhang       }
858cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
859cf1a0672SHong Zhang     }
8601634b032SHong Zhang 
861cf1a0672SHong Zhang     /* set values in C */
862cf1a0672SHong Zhang     row = rstart + i;
863cf1a0672SHong Zhang     api = ptap->api;
864cf1a0672SHong Zhang     apj = ptap->apj;
865cf1a0672SHong Zhang     apJ = apj + api[i];
866cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
867cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
8681634b032SHong Zhang 
869cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
870cf1a0672SHong Zhang     ca = coa + co->i[i];
871cf1a0672SHong Zhang     k  = 0;
872cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++){
873cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
874cf1a0672SHong Zhang       ca[k0]      = apa[apJ[k]];
875cf1a0672SHong Zhang       apa[apJ[k]] = 0.0;
876cf1a0672SHong Zhang       k++;
877cf1a0672SHong Zhang     }
8781634b032SHong Zhang 
879cf1a0672SHong Zhang     /* diagonal part of C */
880cf1a0672SHong Zhang     ca = cda + cd->i[i];
881cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++){
882cf1a0672SHong Zhang       ca[k1]      = apa[apJ[k]];
883cf1a0672SHong Zhang       apa[apJ[k]] = 0.0;
884cf1a0672SHong Zhang       k++;
885cf1a0672SHong Zhang     }
886cf1a0672SHong Zhang 
887cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
888cf1a0672SHong Zhang     ca = coa + co->i[i];
889cf1a0672SHong Zhang     for (; k0<conz; k0++){
890cf1a0672SHong Zhang       ca[k0]      = apa[apJ[k]];
891cf1a0672SHong Zhang       apa[apJ[k]] = 0.0;
892cf1a0672SHong Zhang       k++;
893cf1a0672SHong Zhang     }
894cf1a0672SHong Zhang   }
895cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
896cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
897cf1a0672SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
898cf1a0672SHong Zhang   PetscFunctionReturn(0);
899cf1a0672SHong Zhang }
900cf1a0672SHong Zhang 
901cf1a0672SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */
902cf1a0672SHong Zhang #undef __FUNCT__
903cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
904cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C)
905cf1a0672SHong Zhang {
906cf1a0672SHong Zhang   PetscErrorCode       ierr;
907cf1a0672SHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
908cf1a0672SHong Zhang   Mat                  Cmpi;
909cf1a0672SHong Zhang   Mat_PtAPMPI          *ptap;
910cf1a0672SHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
911cf1a0672SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
912cf1a0672SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
913cf1a0672SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
914cf1a0672SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
915cf1a0672SHong Zhang   PetscInt             i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*nlnk,*lnk,apnz_max=0;
916cf1a0672SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
917cf1a0672SHong Zhang   PetscInt             nlnk_max,armax,prmax;
918cf1a0672SHong Zhang   PetscBT              lnkbt;
919cf1a0672SHong Zhang   PetscReal            afill;
920cf1a0672SHong Zhang   PetscScalar          *apa;
921cf1a0672SHong Zhang   PetscBool            scalable=PETSC_FALSE;
922cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
923cf1a0672SHong Zhang   PetscMPIInt          rank=a->rank;
924cf1a0672SHong Zhang #endif
925cf1a0672SHong Zhang 
926cf1a0672SHong Zhang   PetscFunctionBegin;
927cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
928cf1a0672SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
929cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank);
930cf1a0672SHong Zhang #endif
931cf1a0672SHong Zhang 
932cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
933cf1a0672SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
934cf1a0672SHong Zhang   ptap->apa = PETSC_NULL;
935cf1a0672SHong Zhang 
936cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
937cf1a0672SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
938cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
939cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
940cf1a0672SHong Zhang #endif
941cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
942cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
943cf1a0672SHong Zhang 
944cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
945cf1a0672SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
946cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
947cf1a0672SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
948cf1a0672SHong Zhang 
949cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
950cf1a0672SHong 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);
951cf1a0672SHong Zhang #endif
952cf1a0672SHong Zhang 
953cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
954cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
955cf1a0672SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
956cf1a0672SHong Zhang   ptap->api = api;
957cf1a0672SHong Zhang   api[0]    = 0;
958cf1a0672SHong Zhang 
959cf1a0672SHong Zhang   /* create and initialize a linked list */
960cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
961cf1a0672SHong Zhang   prmax = PetscMax(p_loc->rmax,p_oth->rmax);
962cf1a0672SHong Zhang   nlnk_max = armax*prmax;
963cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
964cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
965cf1a0672SHong 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);
966cf1a0672SHong Zhang #endif
967cf1a0672SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,lnk,nlnk,lnkbt);
968cf1a0672SHong Zhang 
969cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
970cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
971cf1a0672SHong Zhang   current_space = free_space;
972cf1a0672SHong Zhang 
973cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
974cf1a0672SHong Zhang   for (i=0; i<am; i++) {
975cf1a0672SHong Zhang     apnz = 0;
976cf1a0672SHong Zhang     /* diagonal portion of A */
977cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
978cf1a0672SHong Zhang     for (j=0; j<nzi; j++){
979cf1a0672SHong Zhang       row = *adj++;
980cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
981cf1a0672SHong Zhang       Jptr  = pj_loc + pi_loc[row];
982cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
983cf1a0672SHong Zhang       ierr = PetscLLCondensedAddSorted(nlnk_max,pN,pnz,Jptr,nlnk,lnk,lnkbt);CHKERRQ(ierr);
984cf1a0672SHong Zhang     }
985cf1a0672SHong Zhang     /* off-diagonal portion of A */
986cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
987cf1a0672SHong Zhang     for (j=0; j<nzi; j++){
988cf1a0672SHong Zhang       row = *aoj++;
989cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
990cf1a0672SHong Zhang       Jptr  = pj_oth + pi_oth[row];
991cf1a0672SHong Zhang       ierr = PetscLLCondensedAddSorted(nlnk_max,pN,pnz,Jptr,nlnk,lnk,lnkbt);CHKERRQ(ierr);
992cf1a0672SHong Zhang     }
993cf1a0672SHong Zhang 
994cf1a0672SHong Zhang     apnz     = *nlnk;
995cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
996cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
997cf1a0672SHong Zhang 
998cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
999cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
1000cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1001cf1a0672SHong Zhang       nspacedouble++;
1002cf1a0672SHong Zhang     }
1003cf1a0672SHong Zhang 
1004cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
1005cf1a0672SHong Zhang     ierr = PetscLLCondensedClean(nlnk_max,pN,apnz,current_space->array,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1006cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1007cf1a0672SHong Zhang     current_space->array           += apnz;
1008cf1a0672SHong Zhang     current_space->local_used      += apnz;
1009cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
1010cf1a0672SHong Zhang   }
1011cf1a0672SHong Zhang 
1012cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
1013cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
1014cf1a0672SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
1015cf1a0672SHong Zhang   apj  = ptap->apj;
1016cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
1017cf1a0672SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
1018cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
1019cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done..., apnz_max %d\n",rank,apnz_max);
1020cf1a0672SHong Zhang #endif
1021cf1a0672SHong Zhang 
1022cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
1023cf1a0672SHong Zhang   /*----------------------------------------------------*/
1024cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1025cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1026cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1027cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1028cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1029cf1a0672SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1030cf1a0672SHong Zhang 
1031cf1a0672SHong Zhang   /* malloc apa for assembly Cmpi - do I need to do this here? */
1032cf1a0672SHong Zhang   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
1033cf1a0672SHong Zhang   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
1034cf1a0672SHong Zhang   for (i=0; i<am; i++){
1035cf1a0672SHong Zhang     row  = i + rstart;
1036cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
1037cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
1038cf1a0672SHong Zhang     apj += apnz;
1039cf1a0672SHong Zhang   }
1040cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1041cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1042cf1a0672SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
1043cf1a0672SHong Zhang 
1044cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
1045cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
1046cf1a0672SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1047cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
1048cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
1049cf1a0672SHong Zhang 
1050cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1051cf1a0672SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
1052cf1a0672SHong Zhang   c->ptap  = ptap;
1053cf1a0672SHong Zhang 
1054cf1a0672SHong Zhang   *C = Cmpi;
1055cf1a0672SHong Zhang 
1056cf1a0672SHong Zhang   /* set MatInfo */
1057cf1a0672SHong Zhang   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
1058cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
1059cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
1060cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
1061cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
1062cf1a0672SHong Zhang 
1063cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
1064cf1a0672SHong Zhang   if (api[am]) {
1065cf1a0672SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1066cf1a0672SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
1067cf1a0672SHong Zhang   } else {
1068cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1069cf1a0672SHong Zhang   }
1070cf1a0672SHong Zhang #endif
1071cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT)
1072cf1a0672SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank);
1073cf1a0672SHong Zhang #endif
10741634b032SHong Zhang   PetscFunctionReturn(0);
10751634b032SHong Zhang }
1076