1be1d678aSKris Buschelman 22499ec78SHong Zhang /* 32499ec78SHong Zhang Defines matrix-matrix product routines for pairs of MPIAIJ matrices 42499ec78SHong Zhang C = A * B 52499ec78SHong Zhang */ 6*c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 7*c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 8*c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 9*c6db04a5SJed Brown #include <petscbt.h> 10*c6db04a5SJed 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); 39f33d1a9aSHong Zhang if (mult->isrowa){ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);} 40f33d1a9aSHong Zhang if (mult->isrowb){ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);} 41f33d1a9aSHong Zhang if (mult->iscolb){ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);} 42f33d1a9aSHong Zhang if (mult->C_seq){ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);} 43f33d1a9aSHong Zhang if (mult->A_loc){ierr = MatDestroy(mult->A_loc);CHKERRQ(ierr); } 44f33d1a9aSHong Zhang if (mult->B_seq){ierr = MatDestroy(mult->B_seq);CHKERRQ(ierr);} 45f33d1a9aSHong Zhang if (mult->B_loc){ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr);} 46f33d1a9aSHong Zhang if (mult->B_oth){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); 649649163dSHong Zhang if (container) { 65776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 66e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 67dce485f0SBarry Smith A->ops->destroy = mult->destroy; 689649163dSHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 69992edd73SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 70776b82aeSLisandro Dalcin ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 712499ec78SHong Zhang PetscFunctionReturn(0); 722499ec78SHong Zhang } 732499ec78SHong Zhang 742499ec78SHong Zhang #undef __FUNCT__ 75abf897f8SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 76b9c92825SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 77b9c92825SBarry Smith { 78abf897f8SHong Zhang PetscErrorCode ierr; 79abf897f8SHong Zhang Mat_MatMatMultMPI *mult; 80776b82aeSLisandro Dalcin PetscContainer container; 81abf897f8SHong Zhang 82abf897f8SHong Zhang PetscFunctionBegin; 83abf897f8SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 84abf897f8SHong Zhang if (container) { 85776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 86e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 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); 1182499ec78SHong Zhang ierr = MatGetLocalMatCondensed(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); 130f33d1a9aSHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 131776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 132dce485f0SBarry Smith mult->destroy = (*C)->ops->destroy; 133dce485f0SBarry Smith mult->duplicate = (*C)->ops->duplicate; 1342499ec78SHong Zhang 1352499ec78SHong Zhang (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 136abf897f8SHong Zhang (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 1372499ec78SHong Zhang PetscFunctionReturn(0); 1382499ec78SHong Zhang } 1392499ec78SHong Zhang 1402499ec78SHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 1412499ec78SHong Zhang #undef __FUNCT__ 1422499ec78SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 1432499ec78SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 1442499ec78SHong Zhang { 1452499ec78SHong Zhang PetscErrorCode ierr; 1462499ec78SHong Zhang Mat *seq; 1472499ec78SHong Zhang Mat_MatMatMultMPI *mult; 148776b82aeSLisandro Dalcin PetscContainer container; 1492499ec78SHong Zhang 1502499ec78SHong Zhang PetscFunctionBegin; 151f33d1a9aSHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 1522499ec78SHong Zhang if (container) { 153776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 154e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 1552499ec78SHong Zhang 1562499ec78SHong Zhang seq = &mult->B_seq; 1572499ec78SHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1582499ec78SHong Zhang mult->B_seq = *seq; 1592499ec78SHong Zhang 1602499ec78SHong Zhang seq = &mult->A_loc; 1612499ec78SHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1622499ec78SHong Zhang mult->A_loc = *seq; 1632499ec78SHong Zhang 1642499ec78SHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 1652499ec78SHong Zhang 1662499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 167d0f46423SBarry Smith ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 1682499ec78SHong Zhang PetscFunctionReturn(0); 1692499ec78SHong Zhang } 1709123193aSHong Zhang 1719123193aSHong Zhang #undef __FUNCT__ 1729123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 1739123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1749123193aSHong Zhang { 1759123193aSHong Zhang PetscErrorCode ierr; 1769123193aSHong Zhang 1779123193aSHong Zhang PetscFunctionBegin; 1789123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 1799123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 1809123193aSHong Zhang } 1819123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 1829123193aSHong Zhang PetscFunctionReturn(0); 1839123193aSHong Zhang } 1849123193aSHong Zhang 1854324174eSBarry Smith typedef struct { 1864324174eSBarry Smith Mat workB; 1874324174eSBarry Smith PetscScalar *rvalues,*svalues; 1884324174eSBarry Smith MPI_Request *rwaits,*swaits; 1894324174eSBarry Smith } MPIAIJ_MPIDense; 1904324174eSBarry Smith 1917af9e4fdSHong Zhang #undef __FUNCT__ 1927af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 1934324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 1944324174eSBarry Smith { 1954324174eSBarry Smith MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 1964324174eSBarry Smith PetscErrorCode ierr; 1974324174eSBarry Smith 1984324174eSBarry Smith PetscFunctionBegin; 1994324174eSBarry Smith if (contents->workB) {ierr = MatDestroy(contents->workB);CHKERRQ(ierr);} 2004324174eSBarry Smith ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 2014324174eSBarry Smith ierr = PetscFree(contents);CHKERRQ(ierr); 2024324174eSBarry Smith PetscFunctionReturn(0); 2034324174eSBarry Smith } 2044324174eSBarry Smith 2059123193aSHong Zhang #undef __FUNCT__ 2069123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 2079123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 2089123193aSHong Zhang { 2099123193aSHong Zhang PetscErrorCode ierr; 210f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 211d0f46423SBarry Smith PetscInt nz = aij->B->cmap->n; 212776b82aeSLisandro Dalcin PetscContainer cont; 2134324174eSBarry Smith MPIAIJ_MPIDense *contents; 2144324174eSBarry Smith VecScatter ctx = aij->Mvctx; 2154324174eSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 2164324174eSBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 217d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 2189123193aSHong Zhang 2199123193aSHong Zhang PetscFunctionBegin; 220cb2480eaSBarry Smith ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 221d0f46423SBarry Smith ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 222cb2480eaSBarry Smith ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 223cb2480eaSBarry Smith ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 224cb2480eaSBarry Smith ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 225f32d5d43SBarry Smith 2267adad957SLisandro Dalcin ierr = PetscContainerCreate(((PetscObject)A)->comm,&cont);CHKERRQ(ierr); 2274324174eSBarry Smith ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 228776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(cont,contents);CHKERRQ(ierr); 229776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(cont,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 2304324174eSBarry Smith 231f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 232d0f46423SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 2334324174eSBarry Smith 2344324174eSBarry Smith /* Create work arrays needed */ 235d0f46423SBarry Smith ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 236d0f46423SBarry Smith B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 2374324174eSBarry Smith from->n,MPI_Request,&contents->rwaits, 2384324174eSBarry Smith to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 2394324174eSBarry Smith 2404324174eSBarry Smith ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)cont);CHKERRQ(ierr); 241776b82aeSLisandro Dalcin ierr = PetscContainerDestroy(cont);CHKERRQ(ierr); 2429123193aSHong Zhang PetscFunctionReturn(0); 2439123193aSHong Zhang } 2449123193aSHong Zhang 2457af9e4fdSHong Zhang #undef __FUNCT__ 2467af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter" 247f32d5d43SBarry Smith /* 2482636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 2492636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 250f32d5d43SBarry Smith */ 2514324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 252f32d5d43SBarry Smith { 253f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 254f32d5d43SBarry Smith PetscErrorCode ierr; 255f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 256f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 257f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 258f32d5d43SBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 259f32d5d43SBarry Smith PetscInt i,j,k; 260aa5bb8c0SSatish Balay PetscInt *sindices,*sstarts,*rindices,*rstarts; 261aa5bb8c0SSatish Balay PetscMPIInt *sprocs,*rprocs,nrecvs; 262f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 2637adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)A)->comm; 264d0f46423SBarry Smith PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 265f32d5d43SBarry Smith MPI_Status status; 2664324174eSBarry Smith MPIAIJ_MPIDense *contents; 267776b82aeSLisandro Dalcin PetscContainer cont; 2684324174eSBarry Smith Mat workB; 269f32d5d43SBarry Smith 270f32d5d43SBarry Smith PetscFunctionBegin; 2714324174eSBarry Smith ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&cont);CHKERRQ(ierr); 272776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(cont,(void**)&contents);CHKERRQ(ierr); 2734324174eSBarry Smith 2744324174eSBarry Smith workB = *outworkB = contents->workB; 275e32f2f54SBarry 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); 276f32d5d43SBarry Smith sindices = to->indices; 277f32d5d43SBarry Smith sstarts = to->starts; 278f32d5d43SBarry Smith sprocs = to->procs; 2794324174eSBarry Smith swaits = contents->swaits; 2804324174eSBarry Smith svalues = contents->svalues; 281f32d5d43SBarry Smith 282f32d5d43SBarry Smith rindices = from->indices; 283f32d5d43SBarry Smith rstarts = from->starts; 284f32d5d43SBarry Smith rprocs = from->procs; 2854324174eSBarry Smith rwaits = contents->rwaits; 2864324174eSBarry Smith rvalues = contents->rvalues; 287f32d5d43SBarry Smith 288f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 289f32d5d43SBarry Smith ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 290f32d5d43SBarry Smith 291f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 292f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 293f32d5d43SBarry Smith } 294f32d5d43SBarry Smith 295f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 296f32d5d43SBarry Smith /* pack a message at a time */ 297f32d5d43SBarry Smith CHKMEMQ; 298f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 299f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 300f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 3012636ff26SBarry Smith } 3022636ff26SBarry Smith } 303f32d5d43SBarry Smith CHKMEMQ; 304f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 305f32d5d43SBarry Smith } 3062636ff26SBarry Smith 307f32d5d43SBarry Smith nrecvs = from->n; 308f32d5d43SBarry Smith while (nrecvs) { 309f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 310f32d5d43SBarry Smith nrecvs--; 311f32d5d43SBarry Smith /* unpack a message at a time */ 312f32d5d43SBarry Smith CHKMEMQ; 313f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 314f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 315f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 3162636ff26SBarry Smith } 3172636ff26SBarry Smith } 318f32d5d43SBarry Smith CHKMEMQ; 319f32d5d43SBarry Smith } 320cb9801acSJed Brown if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 321f32d5d43SBarry Smith 322f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 323f32d5d43SBarry Smith ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 324f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 325f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 326f32d5d43SBarry Smith PetscFunctionReturn(0); 327f32d5d43SBarry Smith } 328f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 329f32d5d43SBarry Smith 3309123193aSHong Zhang #undef __FUNCT__ 3319123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 3329123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 3339123193aSHong Zhang { 3349123193aSHong Zhang PetscErrorCode ierr; 335f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 336f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 337f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 338f32d5d43SBarry Smith Mat workB; 3399123193aSHong Zhang 3409123193aSHong Zhang PetscFunctionBegin; 3419123193aSHong Zhang 342f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 343f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 344f32d5d43SBarry Smith 345f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 3464324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 347f32d5d43SBarry Smith 348f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 349f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 3509123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3519123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3529123193aSHong Zhang PetscFunctionReturn(0); 3539123193aSHong Zhang } 3549123193aSHong Zhang 355