xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
32499ec78SHong Zhang /*
42499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
52499ec78SHong Zhang           C = A * B
62499ec78SHong Zhang */
7*7c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
8*7c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
9*7c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h"
102499ec78SHong Zhang #include "petscbt.h"
11*7c4f633dSBarry Smith #include "../src/mat/impls/dense/mpi/mpidense.h"
122499ec78SHong Zhang 
132499ec78SHong Zhang #undef __FUNCT__
142499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
152499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
162499ec78SHong Zhang {
172499ec78SHong Zhang   PetscErrorCode ierr;
182499ec78SHong Zhang 
192499ec78SHong Zhang   PetscFunctionBegin;
202499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
212499ec78SHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
222499ec78SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
232499ec78SHong Zhang     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
242499ec78SHong Zhang   } else {
252499ec78SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
262499ec78SHong Zhang   }
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;
3805b42c5fSBarry Smith   ierr = PetscFree(mult->startsj);CHKERRQ(ierr);
3905b42c5fSBarry Smith   ierr = PetscFree(mult->bufa);CHKERRQ(ierr);
40f33d1a9aSHong Zhang   if (mult->isrowa){ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);}
41f33d1a9aSHong Zhang   if (mult->isrowb){ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);}
42f33d1a9aSHong Zhang   if (mult->iscolb){ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);}
43f33d1a9aSHong Zhang   if (mult->C_seq){ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);}
44f33d1a9aSHong Zhang   if (mult->A_loc){ierr = MatDestroy(mult->A_loc);CHKERRQ(ierr); }
45f33d1a9aSHong Zhang   if (mult->B_seq){ierr = MatDestroy(mult->B_seq);CHKERRQ(ierr);}
46f33d1a9aSHong Zhang   if (mult->B_loc){ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr);}
47f33d1a9aSHong Zhang   if (mult->B_oth){ierr = MatDestroy(mult->B_oth);CHKERRQ(ierr);}
4805b42c5fSBarry Smith   ierr = PetscFree(mult->abi);CHKERRQ(ierr);
4905b42c5fSBarry Smith   ierr = PetscFree(mult->abj);CHKERRQ(ierr);
50f33d1a9aSHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
51f33d1a9aSHong Zhang   PetscFunctionReturn(0);
52f33d1a9aSHong Zhang }
53f33d1a9aSHong Zhang 
544a7f9380SHong Zhang EXTERN PetscErrorCode MatDestroy_AIJ(Mat);
552499ec78SHong Zhang #undef __FUNCT__
562499ec78SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
572499ec78SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
582499ec78SHong Zhang {
592499ec78SHong Zhang   PetscErrorCode     ierr;
60776b82aeSLisandro Dalcin   PetscContainer     container;
619649163dSHong Zhang   Mat_MatMatMultMPI  *mult=PETSC_NULL;
622499ec78SHong Zhang 
632499ec78SHong Zhang   PetscFunctionBegin;
649649163dSHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
659649163dSHong Zhang   if (container) {
66776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
679649163dSHong Zhang   } else {
6816bf3c38SBarry Smith     SETERRQ(PETSC_ERR_PLIB,"Container does not exit");
699649163dSHong Zhang   }
7016bf3c38SBarry Smith   A->ops->destroy = mult->MatDestroy;
719649163dSHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
72992edd73SBarry Smith   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
73776b82aeSLisandro Dalcin   ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
742499ec78SHong Zhang   PetscFunctionReturn(0);
752499ec78SHong Zhang }
762499ec78SHong Zhang 
772499ec78SHong Zhang #undef __FUNCT__
78abf897f8SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
79abf897f8SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) {
80abf897f8SHong Zhang   PetscErrorCode     ierr;
81abf897f8SHong Zhang   Mat_MatMatMultMPI  *mult;
82776b82aeSLisandro Dalcin   PetscContainer     container;
83abf897f8SHong Zhang 
84abf897f8SHong Zhang   PetscFunctionBegin;
85abf897f8SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
86abf897f8SHong Zhang   if (container) {
87776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
88abf897f8SHong Zhang   } else {
89abf897f8SHong Zhang     SETERRQ(PETSC_ERR_PLIB,"Container does not exit");
90abf897f8SHong Zhang   }
915c088420SHong Zhang   /* Note: the container is not duplicated, because it requires deep copying of
925c088420SHong Zhang      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
935c088420SHong Zhang      These data sets are only used for repeated calling of MatMatMultNumeric().
945c088420SHong Zhang      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
95abf897f8SHong Zhang   ierr = (*mult->MatDuplicate)(A,op,M);CHKERRQ(ierr);
96abf897f8SHong Zhang   (*M)->ops->destroy   = mult->MatDestroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
97abf897f8SHong Zhang   (*M)->ops->duplicate = mult->MatDuplicate; /* = MatDuplicate_MPIAIJ */
98abf897f8SHong Zhang   PetscFunctionReturn(0);
99abf897f8SHong Zhang }
100abf897f8SHong Zhang 
101abf897f8SHong Zhang #undef __FUNCT__
1022499ec78SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
1032499ec78SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
1042499ec78SHong Zhang {
1052499ec78SHong Zhang   PetscErrorCode     ierr;
1062499ec78SHong Zhang   PetscInt           start,end;
1072499ec78SHong Zhang   Mat_MatMatMultMPI  *mult;
108776b82aeSLisandro Dalcin   PetscContainer     container;
1092499ec78SHong Zhang 
1102499ec78SHong Zhang   PetscFunctionBegin;
111d0f46423SBarry Smith   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
112d0f46423SBarry Smith     SETERRQ4(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);
1132499ec78SHong Zhang   }
1142499ec78SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
1152499ec78SHong Zhang 
1162499ec78SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
1172499ec78SHong Zhang   ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr);
1182499ec78SHong Zhang 
1192499ec78SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
120d0f46423SBarry Smith   start = A->rmap->rstart; end = A->rmap->rend;
1212499ec78SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
1222499ec78SHong Zhang   ierr = MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr);
1232499ec78SHong Zhang 
1242499ec78SHong Zhang   /* compute C_seq = A_seq * B_seq */
1252499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
1262499ec78SHong Zhang 
1272499ec78SHong Zhang   /* create mpi matrix C by concatinating C_seq */
1282499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
129d0f46423SBarry Smith   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
1302499ec78SHong Zhang 
1312499ec78SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
132776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
133776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
134f33d1a9aSHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
135776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
1364a7f9380SHong Zhang   mult->MatDestroy   = (*C)->ops->destroy;
137abf897f8SHong Zhang   mult->MatDuplicate = (*C)->ops->duplicate;
1382499ec78SHong Zhang 
1392499ec78SHong Zhang   (*C)->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
140abf897f8SHong Zhang   (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
1412499ec78SHong Zhang   PetscFunctionReturn(0);
1422499ec78SHong Zhang }
1432499ec78SHong Zhang 
1442499ec78SHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
1452499ec78SHong Zhang #undef __FUNCT__
1462499ec78SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
1472499ec78SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
1482499ec78SHong Zhang {
1492499ec78SHong Zhang   PetscErrorCode     ierr;
1502499ec78SHong Zhang   Mat                *seq;
1512499ec78SHong Zhang   Mat_MatMatMultMPI  *mult;
152776b82aeSLisandro Dalcin   PetscContainer     container;
1532499ec78SHong Zhang 
1542499ec78SHong Zhang   PetscFunctionBegin;
155f33d1a9aSHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
1562499ec78SHong Zhang   if (container) {
157776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
15816bf3c38SBarry Smith   } else {
15916bf3c38SBarry Smith     SETERRQ(PETSC_ERR_PLIB,"Container does not exit");
1602499ec78SHong Zhang   }
1612499ec78SHong Zhang 
1622499ec78SHong Zhang   seq = &mult->B_seq;
1632499ec78SHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
1642499ec78SHong Zhang   mult->B_seq = *seq;
1652499ec78SHong Zhang 
1662499ec78SHong Zhang   seq = &mult->A_loc;
1672499ec78SHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
1682499ec78SHong Zhang   mult->A_loc = *seq;
1692499ec78SHong Zhang 
1702499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
1712499ec78SHong Zhang 
1722499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
173d0f46423SBarry Smith   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
1742499ec78SHong Zhang   PetscFunctionReturn(0);
1752499ec78SHong Zhang }
1769123193aSHong Zhang 
1779123193aSHong Zhang #undef __FUNCT__
1789123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
1799123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1809123193aSHong Zhang {
1819123193aSHong Zhang   PetscErrorCode ierr;
1829123193aSHong Zhang 
1839123193aSHong Zhang   PetscFunctionBegin;
1849123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
1859123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
1869123193aSHong Zhang   }
1879123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
1889123193aSHong Zhang   PetscFunctionReturn(0);
1899123193aSHong Zhang }
1909123193aSHong Zhang 
1914324174eSBarry Smith typedef struct {
1924324174eSBarry Smith   Mat         workB;
1934324174eSBarry Smith   PetscScalar *rvalues,*svalues;
1944324174eSBarry Smith   MPI_Request *rwaits,*swaits;
1954324174eSBarry Smith } MPIAIJ_MPIDense;
1964324174eSBarry Smith 
1974324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
1984324174eSBarry Smith {
1994324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
2004324174eSBarry Smith   PetscErrorCode  ierr;
2014324174eSBarry Smith 
2024324174eSBarry Smith   PetscFunctionBegin;
2034324174eSBarry Smith   if (contents->workB) {ierr = MatDestroy(contents->workB);CHKERRQ(ierr);}
2044324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
2054324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
2064324174eSBarry Smith   PetscFunctionReturn(0);
2074324174eSBarry Smith }
2084324174eSBarry Smith 
2099123193aSHong Zhang #undef __FUNCT__
2109123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
2119123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
2129123193aSHong Zhang {
2139123193aSHong Zhang   PetscErrorCode         ierr;
214f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
215d0f46423SBarry Smith   PetscInt               nz = aij->B->cmap->n;
216776b82aeSLisandro Dalcin   PetscContainer         cont;
2174324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
2184324174eSBarry Smith   VecScatter             ctx = aij->Mvctx;
2194324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
2204324174eSBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
221d0f46423SBarry Smith   PetscInt               m=A->rmap->n,n=B->cmap->n;
2229123193aSHong Zhang 
2239123193aSHong Zhang   PetscFunctionBegin;
224cb2480eaSBarry Smith   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
225d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
226cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
227cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
228cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
229f32d5d43SBarry Smith 
2307adad957SLisandro Dalcin   ierr = PetscContainerCreate(((PetscObject)A)->comm,&cont);CHKERRQ(ierr);
2314324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
232776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(cont,contents);CHKERRQ(ierr);
233776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(cont,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
2344324174eSBarry Smith 
235f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
236d0f46423SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
2374324174eSBarry Smith 
2384324174eSBarry Smith   /* Create work arrays needed */
239d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
240d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
2414324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
2424324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
2434324174eSBarry Smith 
2444324174eSBarry Smith   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)cont);CHKERRQ(ierr);
245776b82aeSLisandro Dalcin   ierr = PetscContainerDestroy(cont);CHKERRQ(ierr);
2469123193aSHong Zhang   PetscFunctionReturn(0);
2479123193aSHong Zhang }
2489123193aSHong Zhang 
249f32d5d43SBarry Smith /*
2502636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
2512636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
252f32d5d43SBarry Smith */
2534324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
254f32d5d43SBarry Smith {
255f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
256f32d5d43SBarry Smith   PetscErrorCode         ierr;
257f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
258f32d5d43SBarry Smith   VecScatter             ctx = aij->Mvctx;
259f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
260f32d5d43SBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
261f32d5d43SBarry Smith   PetscInt               i,j,k;
262aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
263aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
264f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
2657adad957SLisandro Dalcin   MPI_Comm               comm = ((PetscObject)A)->comm;
266d0f46423SBarry Smith   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
267f32d5d43SBarry Smith   MPI_Status             status;
2684324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
269776b82aeSLisandro Dalcin   PetscContainer         cont;
2704324174eSBarry Smith   Mat                    workB;
271f32d5d43SBarry Smith 
272f32d5d43SBarry Smith   PetscFunctionBegin;
2734324174eSBarry Smith   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&cont);CHKERRQ(ierr);
274776b82aeSLisandro Dalcin   ierr = PetscContainerGetPointer(cont,(void**)&contents);CHKERRQ(ierr);
2754324174eSBarry Smith 
2764324174eSBarry Smith   workB = *outworkB = contents->workB;
277d0f46423SBarry Smith   if (nrows != workB->rmap->n) SETERRQ2(PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
278f32d5d43SBarry Smith   sindices  = to->indices;
279f32d5d43SBarry Smith   sstarts   = to->starts;
280f32d5d43SBarry Smith   sprocs    = to->procs;
2814324174eSBarry Smith   swaits    = contents->swaits;
2824324174eSBarry Smith   svalues   = contents->svalues;
283f32d5d43SBarry Smith 
284f32d5d43SBarry Smith   rindices  = from->indices;
285f32d5d43SBarry Smith   rstarts   = from->starts;
286f32d5d43SBarry Smith   rprocs    = from->procs;
2874324174eSBarry Smith   rwaits    = contents->rwaits;
2884324174eSBarry Smith   rvalues   = contents->rvalues;
289f32d5d43SBarry Smith 
290f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
291f32d5d43SBarry Smith   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
292f32d5d43SBarry Smith 
293f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
294f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
295f32d5d43SBarry Smith   }
296f32d5d43SBarry Smith 
297f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
298f32d5d43SBarry Smith     /* pack a message at a time */
299f32d5d43SBarry Smith     CHKMEMQ;
300f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
301f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
302f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
3032636ff26SBarry Smith       }
3042636ff26SBarry Smith     }
305f32d5d43SBarry Smith     CHKMEMQ;
306f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
307f32d5d43SBarry Smith   }
3082636ff26SBarry Smith 
309f32d5d43SBarry Smith   nrecvs = from->n;
310f32d5d43SBarry Smith   while (nrecvs) {
311f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
312f32d5d43SBarry Smith     nrecvs--;
313f32d5d43SBarry Smith     /* unpack a message at a time */
314f32d5d43SBarry Smith     CHKMEMQ;
315f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
316f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
317f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
3182636ff26SBarry Smith       }
3192636ff26SBarry Smith     }
320f32d5d43SBarry Smith     CHKMEMQ;
321f32d5d43SBarry Smith   }
322f32d5d43SBarry Smith   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr)}
323f32d5d43SBarry Smith 
324f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
325f32d5d43SBarry Smith   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
326f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
327f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
328f32d5d43SBarry Smith   PetscFunctionReturn(0);
329f32d5d43SBarry Smith }
330f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
331f32d5d43SBarry Smith 
3329123193aSHong Zhang #undef __FUNCT__
3339123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
3349123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
3359123193aSHong Zhang {
3369123193aSHong Zhang   PetscErrorCode       ierr;
337f32d5d43SBarry Smith   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
338f32d5d43SBarry Smith   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
339f32d5d43SBarry Smith   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
340f32d5d43SBarry Smith   Mat                  workB;
3419123193aSHong Zhang 
3429123193aSHong Zhang   PetscFunctionBegin;
3439123193aSHong Zhang 
344f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
345f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
346f32d5d43SBarry Smith 
347f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
3484324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
349f32d5d43SBarry Smith 
350f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
351f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
3529123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3539123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3549123193aSHong Zhang   PetscFunctionReturn(0);
3559123193aSHong Zhang }
3569123193aSHong Zhang 
357