xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision bfcd1a735ee3385f00a8a7f777aee3ec19d95a27)
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>
112499ec78SHong Zhang 
122499ec78SHong Zhang #undef __FUNCT__
132499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
142499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
152499ec78SHong Zhang {
162499ec78SHong Zhang   PetscErrorCode ierr;
172499ec78SHong Zhang 
182499ec78SHong Zhang   PetscFunctionBegin;
192499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
20598bc09dSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
21a1a4e84aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
22598bc09dSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
232499ec78SHong Zhang   }
24598bc09dSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
25598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
26598bc09dSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
272499ec78SHong Zhang   PetscFunctionReturn(0);
282499ec78SHong Zhang }
292499ec78SHong Zhang 
30f33d1a9aSHong Zhang #undef __FUNCT__
31776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI"
32776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr)
33f33d1a9aSHong Zhang {
34f33d1a9aSHong Zhang   PetscErrorCode       ierr;
359649163dSHong Zhang   Mat_MatMatMultMPI    *mult=(Mat_MatMatMultMPI*)ptr;
36f33d1a9aSHong Zhang 
37f33d1a9aSHong Zhang   PetscFunctionBegin;
386bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr);
396bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr);
406bf464f9SBarry Smith   ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr);
416bf464f9SBarry Smith   ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr);
426bf464f9SBarry Smith   ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr);
436bf464f9SBarry Smith   ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr);
44f33d1a9aSHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
45f33d1a9aSHong Zhang   PetscFunctionReturn(0);
46f33d1a9aSHong Zhang }
47f33d1a9aSHong Zhang 
48598bc09dSHong Zhang #undef __FUNCT__
49a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
50a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
51598bc09dSHong Zhang {
52598bc09dSHong Zhang   PetscErrorCode ierr;
53598bc09dSHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
54598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
55598bc09dSHong Zhang 
56598bc09dSHong Zhang   PetscFunctionBegin;
57598bc09dSHong Zhang   ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
58598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
59a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
60a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
61a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
62a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
63598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
64598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
65598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
66598bc09dSHong Zhang   PetscFunctionReturn(0);
67598bc09dSHong Zhang }
68598bc09dSHong Zhang 
692499ec78SHong Zhang #undef __FUNCT__
70a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32"
71a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A)
722499ec78SHong Zhang {
732499ec78SHong Zhang   PetscErrorCode     ierr;
74776b82aeSLisandro Dalcin   PetscContainer     container;
759649163dSHong Zhang   Mat_MatMatMultMPI  *mult=PETSC_NULL;
762499ec78SHong Zhang 
772499ec78SHong Zhang   PetscFunctionBegin;
789649163dSHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
79bf0cc555SLisandro Dalcin   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
80776b82aeSLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
81dce485f0SBarry Smith   A->ops->destroy   = mult->destroy;
82bf0cc555SLisandro Dalcin   A->ops->duplicate = mult->duplicate;
83bf0cc555SLisandro Dalcin   if (A->ops->destroy) {
84992edd73SBarry Smith     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
85bf0cc555SLisandro Dalcin   }
86bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
872499ec78SHong Zhang   PetscFunctionReturn(0);
882499ec78SHong Zhang }
892499ec78SHong Zhang 
902499ec78SHong Zhang #undef __FUNCT__
91a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32"
92a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M)
93b9c92825SBarry Smith {
94abf897f8SHong Zhang   PetscErrorCode     ierr;
95abf897f8SHong Zhang   Mat_MatMatMultMPI  *mult;
96776b82aeSLisandro Dalcin   PetscContainer     container;
97abf897f8SHong Zhang 
98abf897f8SHong Zhang   PetscFunctionBegin;
99abf897f8SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
100bf0cc555SLisandro Dalcin   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
101776b82aeSLisandro Dalcin   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
1025c088420SHong Zhang   /* Note: the container is not duplicated, because it requires deep copying of
1035c088420SHong Zhang      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
1045c088420SHong Zhang      These data sets are only used for repeated calling of MatMatMultNumeric().
1055c088420SHong Zhang      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
106dce485f0SBarry Smith   ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr);
107dce485f0SBarry Smith   (*M)->ops->destroy   = mult->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
108dce485f0SBarry Smith   (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */
109abf897f8SHong Zhang   PetscFunctionReturn(0);
110abf897f8SHong Zhang }
111abf897f8SHong Zhang 
112abf897f8SHong Zhang #undef __FUNCT__
113a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
114a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
1154ae0847bSHong Zhang {
1164ae0847bSHong Zhang   PetscErrorCode     ierr;
1174ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
1184ae0847bSHong Zhang   Mat_PtAPMPI        *ptap=a->ptap;
1194ae0847bSHong Zhang 
1204ae0847bSHong Zhang   PetscFunctionBegin;
1214ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
1224ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
1234ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
1244ae0847bSHong Zhang   PetscFunctionReturn(0);
1254ae0847bSHong Zhang }
1264ae0847bSHong Zhang 
1274ae0847bSHong Zhang #undef __FUNCT__
128a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
129a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
130598bc09dSHong Zhang {
131598bc09dSHong Zhang   PetscErrorCode     ierr;
1324ae0847bSHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
133598bc09dSHong Zhang   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
134598bc09dSHong Zhang   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
135598bc09dSHong Zhang   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
136598bc09dSHong Zhang   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
137598bc09dSHong Zhang   Mat_SeqAIJ         *p_loc,*p_oth;
138598bc09dSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
139598bc09dSHong Zhang   PetscScalar        *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
140598bc09dSHong Zhang   PetscInt           cm=C->rmap->n,anz,pnz;
141598bc09dSHong Zhang   Mat_PtAPMPI        *ptap=c->ptap;
142598bc09dSHong Zhang   PetscInt           *api,*apj,*apJ,cnz,i,j,k,row;
1434ae0847bSHong Zhang   PetscInt           rstart=C->rmap->rstart,cstart=C->cmap->rstart;
144598bc09dSHong Zhang   PetscInt           cdnz,conz,k0,k1;
145598bc09dSHong Zhang 
146598bc09dSHong Zhang   PetscFunctionBegin;
147a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
148598bc09dSHong Zhang   /*-----------------------------------------------------*/
149a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
150a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
151a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
152598bc09dSHong Zhang 
153598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
154598bc09dSHong Zhang   /*----------------------------------------------------------*/
155598bc09dSHong Zhang   /* get data from symbolic products */
156a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
157a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
158598bc09dSHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
159598bc09dSHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
160598bc09dSHong Zhang 
161598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
162598bc09dSHong Zhang   apa = ptap->apa;
163598bc09dSHong Zhang 
164598bc09dSHong Zhang   for (i=0; i<cm; i++) {
165598bc09dSHong Zhang     /* diagonal portion of A */
166598bc09dSHong Zhang     anz = adi[i+1] - adi[i];
167598bc09dSHong Zhang     adj = ad->j + adi[i];
168598bc09dSHong Zhang     ada = ad->a + adi[i];
169598bc09dSHong Zhang     for (j=0; j<anz; j++) {
170598bc09dSHong Zhang       row = adj[j];
171598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
172598bc09dSHong Zhang       pj  = pj_loc + pi_loc[row];
173598bc09dSHong Zhang       pa  = pa_loc + pi_loc[row];
174598bc09dSHong Zhang 
175598bc09dSHong Zhang       /* perform dense axpy */
176598bc09dSHong Zhang       valtmp = ada[j];
177598bc09dSHong Zhang       for (k=0; k<pnz; k++){
178598bc09dSHong Zhang         apa[pj[k]] += valtmp*pa[k];
179598bc09dSHong Zhang       }
180598bc09dSHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
181598bc09dSHong Zhang     }
182598bc09dSHong Zhang 
183598bc09dSHong Zhang     /* off-diagonal portion of A */
184598bc09dSHong Zhang     anz = aoi[i+1] - aoi[i];
185598bc09dSHong Zhang     aoj = ao->j + aoi[i];
186598bc09dSHong Zhang     aoa = ao->a + aoi[i];
187598bc09dSHong Zhang     for (j=0; j<anz; j++) {
188598bc09dSHong Zhang       row = aoj[j];
189598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
190598bc09dSHong Zhang       pj  = pj_oth + pi_oth[row];
191598bc09dSHong Zhang       pa  = pa_oth + pi_oth[row];
192598bc09dSHong Zhang 
193598bc09dSHong Zhang       /* perform dense axpy */
194598bc09dSHong Zhang       valtmp = aoa[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     /* set values in C */
202598bc09dSHong Zhang     row = rstart + i;
203a1a4e84aSHong Zhang     api = ptap->api;
204a1a4e84aSHong Zhang     apj = ptap->apj;
205598bc09dSHong Zhang     apJ = apj + api[i];
206598bc09dSHong Zhang     cnz = api[i+1] - api[i];
207598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
208598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
209598bc09dSHong Zhang 
210598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
211598bc09dSHong Zhang     ca = coa + co->i[i];
212598bc09dSHong Zhang     k  = 0;
213598bc09dSHong Zhang     for (k0=0; k0<conz; k0++){
214598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
215598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
216598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
217598bc09dSHong Zhang       k++;
218598bc09dSHong Zhang     }
219598bc09dSHong Zhang 
220598bc09dSHong Zhang     /* diagonal part of C */
221598bc09dSHong Zhang     ca = cda + cd->i[i];
222598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++){
223598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
224598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
225598bc09dSHong Zhang       k++;
226598bc09dSHong Zhang     }
227598bc09dSHong Zhang 
228598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
229598bc09dSHong Zhang     ca = coa + co->i[i];
230598bc09dSHong Zhang     for (; k0<conz; k0++){
231598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
232598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
233598bc09dSHong Zhang       k++;
234598bc09dSHong Zhang     }
235598bc09dSHong Zhang   }
236598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
237598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
238598bc09dSHong Zhang   PetscFunctionReturn(0);
239598bc09dSHong Zhang }
240598bc09dSHong Zhang 
241598bc09dSHong Zhang #undef __FUNCT__
242a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
243a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
244598bc09dSHong Zhang {
245598bc09dSHong Zhang   PetscErrorCode       ierr;
246598bc09dSHong Zhang   MPI_Comm             comm=((PetscObject)A)->comm;
247bfcd1a73SHong Zhang   Mat                  Cmpi;
248598bc09dSHong Zhang   Mat_PtAPMPI          *ptap;
249598bc09dSHong Zhang   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
2504ae0847bSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
251bfcd1a73SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2524ae0847bSHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2534ae0847bSHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
254bfcd1a73SHong Zhang   PetscInt             nlnk,*lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
255bfcd1a73SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
256598bc09dSHong Zhang   PetscBT              lnkbt;
257598bc09dSHong Zhang   PetscScalar          *apa;
258bfcd1a73SHong Zhang   PetscReal            afill;
259a1a4e84aSHong Zhang   PetscBool            matmatmult_old=PETSC_FALSE;
260598bc09dSHong Zhang 
261598bc09dSHong Zhang   PetscFunctionBegin;
262a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
263a1a4e84aSHong Zhang     SETERRQ4(PETSC_COMM_SELF,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);
264a1a4e84aSHong Zhang   }
265a1a4e84aSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_old",&matmatmult_old,PETSC_NULL);CHKERRQ(ierr);
266a1a4e84aSHong Zhang   if (matmatmult_old){
267a1a4e84aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_32(A,P,fill,C);;CHKERRQ(ierr);
268a1a4e84aSHong Zhang     PetscFunctionReturn(0);
269a1a4e84aSHong Zhang   }
270a1a4e84aSHong Zhang 
271598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
272598bc09dSHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
273598bc09dSHong Zhang   ptap->abnz_max = 0;
274598bc09dSHong Zhang 
275598bc09dSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
276bfcd1a73SHong Zhang   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
277bfcd1a73SHong Zhang   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
278598bc09dSHong Zhang   ptap->apa = apa;
279598bc09dSHong Zhang 
280598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
281a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
282598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
283a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
284598bc09dSHong Zhang 
285a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
286a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
287598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
288598bc09dSHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
289598bc09dSHong Zhang 
290598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
291598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
292598bc09dSHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
293a1a4e84aSHong Zhang   ptap->api = api;
294598bc09dSHong Zhang   api[0]    = 0;
295598bc09dSHong Zhang 
296598bc09dSHong Zhang   /* create and initialize a linked list */
297598bc09dSHong Zhang   nlnk = pN+1;
298598bc09dSHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
299598bc09dSHong Zhang 
300bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
301bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
302598bc09dSHong Zhang   current_space = free_space;
303598bc09dSHong Zhang 
304598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
305598bc09dSHong Zhang   for (i=0; i<am; i++) {
306598bc09dSHong Zhang     apnz = 0;
307598bc09dSHong Zhang     /* diagonal portion of A */
308598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
309598bc09dSHong Zhang     for (j=0; j<nzi; j++){
310598bc09dSHong Zhang       row = *adj++;
311598bc09dSHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
312598bc09dSHong Zhang       Jptr  = pj_loc + pi_loc[row];
313598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
314598bc09dSHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
315598bc09dSHong Zhang       apnz += nlnk;
316598bc09dSHong Zhang     }
317598bc09dSHong Zhang     /* off-diagonal portion of A */
318598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
319598bc09dSHong Zhang     for (j=0; j<nzi; j++){
320598bc09dSHong Zhang       row = *aoj++;
321598bc09dSHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
322598bc09dSHong Zhang       Jptr  = pj_oth + pi_oth[row];
323598bc09dSHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
324598bc09dSHong Zhang       apnz += nlnk;
325598bc09dSHong Zhang     }
326598bc09dSHong Zhang 
327598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
328598bc09dSHong Zhang     if (ptap->abnz_max < apnz) ptap->abnz_max = apnz;
329598bc09dSHong Zhang 
330598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
331598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
332598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
333598bc09dSHong Zhang       nspacedouble++;
334598bc09dSHong Zhang     }
335598bc09dSHong Zhang 
336598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
337598bc09dSHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
338598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
339598bc09dSHong Zhang     current_space->array           += apnz;
340598bc09dSHong Zhang     current_space->local_used      += apnz;
341598bc09dSHong Zhang     current_space->local_remaining -= apnz;
342598bc09dSHong Zhang   }
343598bc09dSHong Zhang 
344598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
345598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
346a1a4e84aSHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
347a1a4e84aSHong Zhang   apj  = ptap->apj;
348a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
349598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
350598bc09dSHong Zhang 
351bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
352598bc09dSHong Zhang   /*----------------------------------------------------*/
353bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
354bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
355bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
356bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
357598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
358bfcd1a73SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
359598bc09dSHong Zhang   for (i=0; i<am; i++){
360598bc09dSHong Zhang     row  = i + rstart;
361598bc09dSHong Zhang     apnz = api[i+1] - api[i];
362bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
363598bc09dSHong Zhang     apj += apnz;
364598bc09dSHong Zhang   }
365bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
366bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
367598bc09dSHong Zhang 
368bfcd1a73SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
369bfcd1a73SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
370bfcd1a73SHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
371bfcd1a73SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
372bfcd1a73SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
373598bc09dSHong Zhang 
374bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
375bfcd1a73SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
376598bc09dSHong Zhang   c->ptap  = ptap;
377598bc09dSHong Zhang 
378bfcd1a73SHong Zhang   *C = Cmpi;
379bfcd1a73SHong Zhang 
380bfcd1a73SHong Zhang   /* set MatInfo */
381bfcd1a73SHong Zhang   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
382bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
383bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
384bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
385bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
386bfcd1a73SHong Zhang 
387bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
388bfcd1a73SHong Zhang   if (api[am]) {
389bfcd1a73SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
390bfcd1a73SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
391bfcd1a73SHong Zhang   } else {
392bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
393bfcd1a73SHong Zhang   }
394bfcd1a73SHong Zhang #endif
395598bc09dSHong Zhang   PetscFunctionReturn(0);
396598bc09dSHong Zhang }
397598bc09dSHong Zhang 
398a1a4e84aSHong Zhang /* implementation used in PETSc-3.2 */
399a1a4e84aSHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
400598bc09dSHong Zhang #undef __FUNCT__
401a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_32"
402a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_32(Mat A,Mat B,Mat C)
403a1a4e84aSHong Zhang {
404a1a4e84aSHong Zhang   PetscErrorCode     ierr;
405a1a4e84aSHong Zhang   Mat                *seq;
406a1a4e84aSHong Zhang   Mat_MatMatMultMPI  *mult;
407a1a4e84aSHong Zhang   PetscContainer     container;
408a1a4e84aSHong Zhang 
409a1a4e84aSHong Zhang   PetscFunctionBegin;
410a1a4e84aSHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
411a1a4e84aSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
412a1a4e84aSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
413a1a4e84aSHong Zhang   seq = &mult->B_seq;
414a1a4e84aSHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
415a1a4e84aSHong Zhang   mult->B_seq = *seq;
416a1a4e84aSHong Zhang 
417a1a4e84aSHong Zhang   seq = &mult->A_loc;
418a1a4e84aSHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
419a1a4e84aSHong Zhang   mult->A_loc = *seq;
420a1a4e84aSHong Zhang 
421a1a4e84aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
422a1a4e84aSHong Zhang 
423a1a4e84aSHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
424a1a4e84aSHong Zhang   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
425a1a4e84aSHong Zhang   PetscFunctionReturn(0);
426a1a4e84aSHong Zhang }
427a1a4e84aSHong Zhang 
428a1a4e84aSHong Zhang #undef __FUNCT__
429a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_32"
430a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_32(Mat A,Mat B,PetscReal fill,Mat *C)
4312499ec78SHong Zhang {
4322499ec78SHong Zhang   PetscErrorCode     ierr;
4332499ec78SHong Zhang   Mat_MatMatMultMPI  *mult;
434776b82aeSLisandro Dalcin   PetscContainer     container;
435d5821860SHong Zhang   Mat                AB,*seq;
436d5821860SHong Zhang   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
437d5821860SHong Zhang   PetscInt           *idx,i,start,ncols,nzA,nzB,*cmap,imark;
4382499ec78SHong Zhang 
4392499ec78SHong Zhang   PetscFunctionBegin;
440d0f46423SBarry Smith   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
441e32f2f54SBarry Smith     SETERRQ4(PETSC_COMM_SELF,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);
4422499ec78SHong Zhang   }
443598bc09dSHong Zhang 
4442499ec78SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
4452499ec78SHong Zhang 
446d5821860SHong Zhang   /* get isrowb: nonzero col of A */
447d5821860SHong Zhang   start = A->cmap->rstart;
448d5821860SHong Zhang   cmap  = a->garray;
449d5821860SHong Zhang   nzA   = a->A->cmap->n;
450d5821860SHong Zhang   nzB   = a->B->cmap->n;
451d5821860SHong Zhang   ierr  = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr);
452d5821860SHong Zhang   ncols = 0;
453d5821860SHong Zhang   for (i=0; i<nzB; i++) {  /* row < local row index */
454d5821860SHong Zhang     if (cmap[i] < start) idx[ncols++] = cmap[i];
455d5821860SHong Zhang     else break;
456d5821860SHong Zhang   }
457d5821860SHong Zhang   imark = i;
458d5821860SHong Zhang   for (i=0; i<nzA; i++) idx[ncols++] = start + i;  /* local rows */
459d5821860SHong Zhang   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */
460d5821860SHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr);
461d5821860SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr);
462d5821860SHong Zhang 
463d5821860SHong Zhang   /*  get isrowa: all local rows of A */
464d5821860SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr);
465d5821860SHong Zhang 
466d5821860SHong Zhang   /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */
4672499ec78SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
468d5821860SHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
469d5821860SHong Zhang   mult->B_seq = *seq;
470d5821860SHong Zhang   ierr = PetscFree(seq);CHKERRQ(ierr);
4712499ec78SHong Zhang 
4722499ec78SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
473d5821860SHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
474d5821860SHong Zhang   mult->A_loc = *seq;
475d5821860SHong Zhang   ierr = PetscFree(seq);CHKERRQ(ierr);
4762499ec78SHong Zhang 
4772499ec78SHong Zhang   /* compute C_seq = A_seq * B_seq */
478*9b8102ccSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr);
479*9b8102ccSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
4802499ec78SHong Zhang 
4812499ec78SHong Zhang   /* create mpi matrix C by concatinating C_seq */
4822499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
483*9b8102ccSHong Zhang   ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr);
484*9b8102ccSHong Zhang   ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr);
4852499ec78SHong Zhang 
4862499ec78SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
487776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
488776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
489776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
490d5821860SHong Zhang   ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
491b160f555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
492d5821860SHong Zhang   mult->destroy      = AB->ops->destroy;
493d5821860SHong Zhang   mult->duplicate    = AB->ops->duplicate;
494a1a4e84aSHong Zhang   AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_32;
495a1a4e84aSHong Zhang   AB->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult_32;
496a1a4e84aSHong Zhang   AB->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult_32;
4978cdbd757SHong Zhang   AB->ops->matmult        = MatMatMult_MPIAIJ_MPIAIJ;
498d5821860SHong Zhang   *C                      = AB;
4992499ec78SHong Zhang   PetscFunctionReturn(0);
5002499ec78SHong Zhang }
5012499ec78SHong Zhang 
5029123193aSHong Zhang #undef __FUNCT__
5039123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
5049123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
5059123193aSHong Zhang {
5069123193aSHong Zhang   PetscErrorCode ierr;
5079123193aSHong Zhang 
5089123193aSHong Zhang   PetscFunctionBegin;
5099123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
5109123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
5119123193aSHong Zhang   }
5129123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
5139123193aSHong Zhang   PetscFunctionReturn(0);
5149123193aSHong Zhang }
5159123193aSHong Zhang 
5164324174eSBarry Smith typedef struct {
5174324174eSBarry Smith   Mat         workB;
5184324174eSBarry Smith   PetscScalar *rvalues,*svalues;
5194324174eSBarry Smith   MPI_Request *rwaits,*swaits;
5204324174eSBarry Smith } MPIAIJ_MPIDense;
5214324174eSBarry Smith 
5227af9e4fdSHong Zhang #undef __FUNCT__
5237af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
5244324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
5254324174eSBarry Smith {
5264324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
5274324174eSBarry Smith   PetscErrorCode  ierr;
5284324174eSBarry Smith 
5294324174eSBarry Smith   PetscFunctionBegin;
5306bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
5314324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
5324324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
5334324174eSBarry Smith   PetscFunctionReturn(0);
5344324174eSBarry Smith }
5354324174eSBarry Smith 
5369123193aSHong Zhang #undef __FUNCT__
5379123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
5389123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
5399123193aSHong Zhang {
5409123193aSHong Zhang   PetscErrorCode         ierr;
541f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
542d0f46423SBarry Smith   PetscInt               nz = aij->B->cmap->n;
543bf0cc555SLisandro Dalcin   PetscContainer         container;
5444324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
5454324174eSBarry Smith   VecScatter             ctx = aij->Mvctx;
5464324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
5474324174eSBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
548d0f46423SBarry Smith   PetscInt               m=A->rmap->n,n=B->cmap->n;
5499123193aSHong Zhang 
5509123193aSHong Zhang   PetscFunctionBegin;
551cb2480eaSBarry Smith   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
552d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
553cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
554cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
555cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5568cdbd757SHong Zhang   (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense;
557f32d5d43SBarry Smith 
5584324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
559f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
560d0f46423SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
5614324174eSBarry Smith   /* Create work arrays needed */
562d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
563d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
5644324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
5654324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
5664324174eSBarry Smith 
567bf0cc555SLisandro Dalcin   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
568bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
569bf0cc555SLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
570bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
571bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
5729123193aSHong Zhang   PetscFunctionReturn(0);
5739123193aSHong Zhang }
5749123193aSHong Zhang 
5757af9e4fdSHong Zhang #undef __FUNCT__
5767af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
577f32d5d43SBarry Smith /*
5782636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
5792636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
580f32d5d43SBarry Smith */
5814324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
582f32d5d43SBarry Smith {
583f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
584f32d5d43SBarry Smith   PetscErrorCode         ierr;
585f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
586f32d5d43SBarry Smith   VecScatter             ctx = aij->Mvctx;
587f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
588f32d5d43SBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
589f32d5d43SBarry Smith   PetscInt               i,j,k;
590aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
591aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
592f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
5937adad957SLisandro Dalcin   MPI_Comm               comm = ((PetscObject)A)->comm;
594d0f46423SBarry Smith   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
595f32d5d43SBarry Smith   MPI_Status             status;
5964324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
597bf0cc555SLisandro Dalcin   PetscContainer         container;
5984324174eSBarry Smith   Mat                    workB;
599f32d5d43SBarry Smith 
600f32d5d43SBarry Smith   PetscFunctionBegin;
601bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
602bf0cc555SLisandro Dalcin   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
603bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
6044324174eSBarry Smith 
6054324174eSBarry Smith   workB = *outworkB = contents->workB;
606e32f2f54SBarry Smith   if (nrows != workB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
607f32d5d43SBarry Smith   sindices  = to->indices;
608f32d5d43SBarry Smith   sstarts   = to->starts;
609f32d5d43SBarry Smith   sprocs    = to->procs;
6104324174eSBarry Smith   swaits    = contents->swaits;
6114324174eSBarry Smith   svalues   = contents->svalues;
612f32d5d43SBarry Smith 
613f32d5d43SBarry Smith   rindices  = from->indices;
614f32d5d43SBarry Smith   rstarts   = from->starts;
615f32d5d43SBarry Smith   rprocs    = from->procs;
6164324174eSBarry Smith   rwaits    = contents->rwaits;
6174324174eSBarry Smith   rvalues   = contents->rvalues;
618f32d5d43SBarry Smith 
619f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
620f32d5d43SBarry Smith   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
621f32d5d43SBarry Smith 
622f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
623f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
624f32d5d43SBarry Smith   }
625f32d5d43SBarry Smith 
626f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
627f32d5d43SBarry Smith     /* pack a message at a time */
628f32d5d43SBarry Smith     CHKMEMQ;
629f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
630f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
631f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
6322636ff26SBarry Smith       }
6332636ff26SBarry Smith     }
634f32d5d43SBarry Smith     CHKMEMQ;
635f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
636f32d5d43SBarry Smith   }
6372636ff26SBarry Smith 
638f32d5d43SBarry Smith   nrecvs = from->n;
639f32d5d43SBarry Smith   while (nrecvs) {
640f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
641f32d5d43SBarry Smith     nrecvs--;
642f32d5d43SBarry Smith     /* unpack a message at a time */
643f32d5d43SBarry Smith     CHKMEMQ;
644f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
645f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
646f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
6472636ff26SBarry Smith       }
6482636ff26SBarry Smith     }
649f32d5d43SBarry Smith     CHKMEMQ;
650f32d5d43SBarry Smith   }
651cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
652f32d5d43SBarry Smith 
653f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
654f32d5d43SBarry Smith   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
655f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
656f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
657f32d5d43SBarry Smith   PetscFunctionReturn(0);
658f32d5d43SBarry Smith }
659f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
660f32d5d43SBarry Smith 
6619123193aSHong Zhang #undef __FUNCT__
6629123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
6639123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
6649123193aSHong Zhang {
6659123193aSHong Zhang   PetscErrorCode       ierr;
666f32d5d43SBarry Smith   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
667f32d5d43SBarry Smith   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
668f32d5d43SBarry Smith   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
669f32d5d43SBarry Smith   Mat                  workB;
6709123193aSHong Zhang 
6719123193aSHong Zhang   PetscFunctionBegin;
6729123193aSHong Zhang 
673f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
674f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
675f32d5d43SBarry Smith 
676f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
6774324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
678f32d5d43SBarry Smith 
679f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
680f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
6819123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6829123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6839123193aSHong Zhang   PetscFunctionReturn(0);
6849123193aSHong Zhang }
6859123193aSHong Zhang 
686