xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision bf0cc55543cd83e035744be2f77202b216d1436e)
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){
202499ec78SHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
212499ec78SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
222499ec78SHong Zhang     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
232499ec78SHong Zhang   } else {
24e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
252499ec78SHong Zhang   }
262499ec78SHong Zhang   PetscFunctionReturn(0);
272499ec78SHong Zhang }
282499ec78SHong Zhang 
29f33d1a9aSHong Zhang #undef __FUNCT__
30776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI"
31776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr)
32f33d1a9aSHong Zhang {
33f33d1a9aSHong Zhang   PetscErrorCode       ierr;
349649163dSHong Zhang   Mat_MatMatMultMPI    *mult=(Mat_MatMatMultMPI*)ptr;
35f33d1a9aSHong Zhang 
36f33d1a9aSHong Zhang   PetscFunctionBegin;
371d79065fSBarry Smith   ierr = PetscFree2(mult->startsj,mult->startsj_r);CHKERRQ(ierr);
3805b42c5fSBarry Smith   ierr = PetscFree(mult->bufa);CHKERRQ(ierr);
396bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr);
406bf464f9SBarry Smith   ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr);
416bf464f9SBarry Smith   ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr);
426bf464f9SBarry Smith   ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr);
436bf464f9SBarry Smith   ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr);
446bf464f9SBarry Smith   ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr);
456bf464f9SBarry Smith   ierr = MatDestroy(&mult->B_loc);CHKERRQ(ierr);
466bf464f9SBarry Smith   ierr = MatDestroy(&mult->B_oth);CHKERRQ(ierr);
4705b42c5fSBarry Smith   ierr = PetscFree(mult->abi);CHKERRQ(ierr);
4805b42c5fSBarry Smith   ierr = PetscFree(mult->abj);CHKERRQ(ierr);
49f33d1a9aSHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
50f33d1a9aSHong Zhang   PetscFunctionReturn(0);
51f33d1a9aSHong Zhang }
52f33d1a9aSHong Zhang 
5309573ac7SBarry Smith extern PetscErrorCode MatDestroy_AIJ(Mat);
542499ec78SHong Zhang #undef __FUNCT__
552499ec78SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
562499ec78SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
572499ec78SHong Zhang {
582499ec78SHong Zhang   PetscErrorCode     ierr;
59776b82aeSLisandro Dalcin   PetscContainer     container;
609649163dSHong Zhang   Mat_MatMatMultMPI  *mult=PETSC_NULL;
612499ec78SHong Zhang 
622499ec78SHong Zhang   PetscFunctionBegin;
639649163dSHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
64*bf0cc555SLisandro Dalcin   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
65776b82aeSLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
66dce485f0SBarry Smith   A->ops->destroy   = mult->destroy;
67*bf0cc555SLisandro Dalcin   A->ops->duplicate = mult->duplicate;
68*bf0cc555SLisandro Dalcin   if (A->ops->destroy) {
69992edd73SBarry Smith     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
70*bf0cc555SLisandro Dalcin   }
71*bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
722499ec78SHong Zhang   PetscFunctionReturn(0);
732499ec78SHong Zhang }
742499ec78SHong Zhang 
752499ec78SHong Zhang #undef __FUNCT__
76abf897f8SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
77b9c92825SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
78b9c92825SBarry Smith {
79abf897f8SHong Zhang   PetscErrorCode     ierr;
80abf897f8SHong Zhang   Mat_MatMatMultMPI  *mult;
81776b82aeSLisandro Dalcin   PetscContainer     container;
82abf897f8SHong Zhang 
83abf897f8SHong Zhang   PetscFunctionBegin;
84abf897f8SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
85*bf0cc555SLisandro Dalcin   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
86776b82aeSLisandro Dalcin   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
875c088420SHong Zhang   /* Note: the container is not duplicated, because it requires deep copying of
885c088420SHong Zhang      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
895c088420SHong Zhang      These data sets are only used for repeated calling of MatMatMultNumeric().
905c088420SHong Zhang      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
91dce485f0SBarry Smith   ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr);
92dce485f0SBarry Smith   (*M)->ops->destroy   = mult->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
93dce485f0SBarry Smith   (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */
94abf897f8SHong Zhang   PetscFunctionReturn(0);
95abf897f8SHong Zhang }
96abf897f8SHong Zhang 
97abf897f8SHong Zhang #undef __FUNCT__
982499ec78SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
992499ec78SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
1002499ec78SHong Zhang {
1012499ec78SHong Zhang   PetscErrorCode     ierr;
1022499ec78SHong Zhang   PetscInt           start,end;
1032499ec78SHong Zhang   Mat_MatMatMultMPI  *mult;
104776b82aeSLisandro Dalcin   PetscContainer     container;
1052499ec78SHong Zhang 
1062499ec78SHong Zhang   PetscFunctionBegin;
107d0f46423SBarry Smith   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
108e32f2f54SBarry 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);
1092499ec78SHong Zhang   }
1102499ec78SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
1112499ec78SHong Zhang 
1122499ec78SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
1132499ec78SHong Zhang   ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr);
1142499ec78SHong Zhang 
1152499ec78SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
116d0f46423SBarry Smith   start = A->rmap->rstart; end = A->rmap->rend;
1172499ec78SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
1184a2b5492SBarry Smith   ierr = MatMPIAIJGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr);
1192499ec78SHong Zhang 
1202499ec78SHong Zhang   /* compute C_seq = A_seq * B_seq */
1212499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
1222499ec78SHong Zhang 
1232499ec78SHong Zhang   /* create mpi matrix C by concatinating C_seq */
1242499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
125d0f46423SBarry Smith   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
1262499ec78SHong Zhang 
1272499ec78SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
128776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
129776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
130776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
131*bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
132dce485f0SBarry Smith   mult->destroy   = (*C)->ops->destroy;
133dce485f0SBarry Smith   mult->duplicate = (*C)->ops->duplicate;
1342499ec78SHong Zhang   (*C)->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
135abf897f8SHong Zhang   (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
1362499ec78SHong Zhang   PetscFunctionReturn(0);
1372499ec78SHong Zhang }
1382499ec78SHong Zhang 
1392499ec78SHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
1402499ec78SHong Zhang #undef __FUNCT__
1412499ec78SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
1422499ec78SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
1432499ec78SHong Zhang {
1442499ec78SHong Zhang   PetscErrorCode     ierr;
1452499ec78SHong Zhang   Mat                *seq;
1462499ec78SHong Zhang   Mat_MatMatMultMPI  *mult;
147776b82aeSLisandro Dalcin   PetscContainer     container;
1482499ec78SHong Zhang 
1492499ec78SHong Zhang   PetscFunctionBegin;
150f33d1a9aSHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
151*bf0cc555SLisandro Dalcin   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
152776b82aeSLisandro Dalcin   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
1532499ec78SHong Zhang 
1542499ec78SHong Zhang   seq = &mult->B_seq;
1552499ec78SHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
1562499ec78SHong Zhang   mult->B_seq = *seq;
1572499ec78SHong Zhang 
1582499ec78SHong Zhang   seq = &mult->A_loc;
1592499ec78SHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
1602499ec78SHong Zhang   mult->A_loc = *seq;
1612499ec78SHong Zhang 
1622499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
1632499ec78SHong Zhang 
1642499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
165d0f46423SBarry Smith   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
1662499ec78SHong Zhang   PetscFunctionReturn(0);
1672499ec78SHong Zhang }
1689123193aSHong Zhang 
1699123193aSHong Zhang #undef __FUNCT__
1709123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
1719123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1729123193aSHong Zhang {
1739123193aSHong Zhang   PetscErrorCode ierr;
1749123193aSHong Zhang 
1759123193aSHong Zhang   PetscFunctionBegin;
1769123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
1779123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
1789123193aSHong Zhang   }
1799123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
1809123193aSHong Zhang   PetscFunctionReturn(0);
1819123193aSHong Zhang }
1829123193aSHong Zhang 
1834324174eSBarry Smith typedef struct {
1844324174eSBarry Smith   Mat         workB;
1854324174eSBarry Smith   PetscScalar *rvalues,*svalues;
1864324174eSBarry Smith   MPI_Request *rwaits,*swaits;
1874324174eSBarry Smith } MPIAIJ_MPIDense;
1884324174eSBarry Smith 
1897af9e4fdSHong Zhang #undef __FUNCT__
1907af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
1914324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
1924324174eSBarry Smith {
1934324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
1944324174eSBarry Smith   PetscErrorCode  ierr;
1954324174eSBarry Smith 
1964324174eSBarry Smith   PetscFunctionBegin;
1976bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
1984324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
1994324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
2004324174eSBarry Smith   PetscFunctionReturn(0);
2014324174eSBarry Smith }
2024324174eSBarry Smith 
2039123193aSHong Zhang #undef __FUNCT__
2049123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
2059123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
2069123193aSHong Zhang {
2079123193aSHong Zhang   PetscErrorCode         ierr;
208f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
209d0f46423SBarry Smith   PetscInt               nz = aij->B->cmap->n;
210*bf0cc555SLisandro Dalcin   PetscContainer         container;
2114324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
2124324174eSBarry Smith   VecScatter             ctx = aij->Mvctx;
2134324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
2144324174eSBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
215d0f46423SBarry Smith   PetscInt               m=A->rmap->n,n=B->cmap->n;
2169123193aSHong Zhang 
2179123193aSHong Zhang   PetscFunctionBegin;
218cb2480eaSBarry Smith   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
219d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
220cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
221cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
222cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223f32d5d43SBarry Smith 
2244324174eSBarry Smith   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
225f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
226d0f46423SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
2274324174eSBarry Smith   /* Create work arrays needed */
228d0f46423SBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
229d0f46423SBarry Smith                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
2304324174eSBarry Smith                       from->n,MPI_Request,&contents->rwaits,
2314324174eSBarry Smith                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
2324324174eSBarry Smith 
233*bf0cc555SLisandro Dalcin   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
234*bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
235*bf0cc555SLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
236*bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
237*bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
2389123193aSHong Zhang   PetscFunctionReturn(0);
2399123193aSHong Zhang }
2409123193aSHong Zhang 
2417af9e4fdSHong Zhang #undef __FUNCT__
2427af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
243f32d5d43SBarry Smith /*
2442636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
2452636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
246f32d5d43SBarry Smith */
2474324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
248f32d5d43SBarry Smith {
249f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
250f32d5d43SBarry Smith   PetscErrorCode         ierr;
251f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
252f32d5d43SBarry Smith   VecScatter             ctx = aij->Mvctx;
253f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
254f32d5d43SBarry Smith   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
255f32d5d43SBarry Smith   PetscInt               i,j,k;
256aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
257aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
258f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
2597adad957SLisandro Dalcin   MPI_Comm               comm = ((PetscObject)A)->comm;
260d0f46423SBarry Smith   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
261f32d5d43SBarry Smith   MPI_Status             status;
2624324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
263*bf0cc555SLisandro Dalcin   PetscContainer         container;
2644324174eSBarry Smith   Mat                    workB;
265f32d5d43SBarry Smith 
266f32d5d43SBarry Smith   PetscFunctionBegin;
267*bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
268*bf0cc555SLisandro Dalcin   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
269*bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
2704324174eSBarry Smith 
2714324174eSBarry Smith   workB = *outworkB = contents->workB;
272e32f2f54SBarry 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);
273f32d5d43SBarry Smith   sindices  = to->indices;
274f32d5d43SBarry Smith   sstarts   = to->starts;
275f32d5d43SBarry Smith   sprocs    = to->procs;
2764324174eSBarry Smith   swaits    = contents->swaits;
2774324174eSBarry Smith   svalues   = contents->svalues;
278f32d5d43SBarry Smith 
279f32d5d43SBarry Smith   rindices  = from->indices;
280f32d5d43SBarry Smith   rstarts   = from->starts;
281f32d5d43SBarry Smith   rprocs    = from->procs;
2824324174eSBarry Smith   rwaits    = contents->rwaits;
2834324174eSBarry Smith   rvalues   = contents->rvalues;
284f32d5d43SBarry Smith 
285f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
286f32d5d43SBarry Smith   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
287f32d5d43SBarry Smith 
288f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
289f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
290f32d5d43SBarry Smith   }
291f32d5d43SBarry Smith 
292f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
293f32d5d43SBarry Smith     /* pack a message at a time */
294f32d5d43SBarry Smith     CHKMEMQ;
295f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
296f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
297f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
2982636ff26SBarry Smith       }
2992636ff26SBarry Smith     }
300f32d5d43SBarry Smith     CHKMEMQ;
301f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
302f32d5d43SBarry Smith   }
3032636ff26SBarry Smith 
304f32d5d43SBarry Smith   nrecvs = from->n;
305f32d5d43SBarry Smith   while (nrecvs) {
306f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
307f32d5d43SBarry Smith     nrecvs--;
308f32d5d43SBarry Smith     /* unpack a message at a time */
309f32d5d43SBarry Smith     CHKMEMQ;
310f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
311f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
312f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
3132636ff26SBarry Smith       }
3142636ff26SBarry Smith     }
315f32d5d43SBarry Smith     CHKMEMQ;
316f32d5d43SBarry Smith   }
317cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
318f32d5d43SBarry Smith 
319f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
320f32d5d43SBarry Smith   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
321f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
322f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
323f32d5d43SBarry Smith   PetscFunctionReturn(0);
324f32d5d43SBarry Smith }
325f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
326f32d5d43SBarry Smith 
3279123193aSHong Zhang #undef __FUNCT__
3289123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
3299123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
3309123193aSHong Zhang {
3319123193aSHong Zhang   PetscErrorCode       ierr;
332f32d5d43SBarry Smith   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
333f32d5d43SBarry Smith   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
334f32d5d43SBarry Smith   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
335f32d5d43SBarry Smith   Mat                  workB;
3369123193aSHong Zhang 
3379123193aSHong Zhang   PetscFunctionBegin;
3389123193aSHong Zhang 
339f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
340f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
341f32d5d43SBarry Smith 
342f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
3434324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
344f32d5d43SBarry Smith 
345f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
346f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
3479123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3489123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3499123193aSHong Zhang   PetscFunctionReturn(0);
3509123193aSHong Zhang }
3519123193aSHong Zhang 
352