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 */ 72499ec78SHong Zhang #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 82499ec78SHong Zhang #include "src/mat/utils/freespace.h" 92499ec78SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 102499ec78SHong Zhang #include "petscbt.h" 114ae313f4SHong Zhang #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__ 31f33d1a9aSHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_Mat_MatMatMultMPI" 329649163dSHong Zhang PetscErrorCode PetscObjectContainerDestroy_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; 609649163dSHong Zhang PetscObjectContainer 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) { 669649163dSHong Zhang ierr = PetscObjectContainerGetPointer(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); 739649163dSHong Zhang ierr = PetscObjectContainerDestroy(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; 82abf897f8SHong Zhang PetscObjectContainer container; 83abf897f8SHong Zhang 84abf897f8SHong Zhang PetscFunctionBegin; 85abf897f8SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 86abf897f8SHong Zhang if (container) { 87abf897f8SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 88abf897f8SHong Zhang } else { 89abf897f8SHong Zhang SETERRQ(PETSC_ERR_PLIB,"Container does not exit"); 90abf897f8SHong Zhang } 91abf897f8SHong Zhang ierr = (*mult->MatDuplicate)(A,op,M);CHKERRQ(ierr); 92abf897f8SHong Zhang (*M)->ops->destroy = mult->MatDestroy; /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 93abf897f8SHong Zhang (*M)->ops->duplicate = mult->MatDuplicate; /* =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; 1042499ec78SHong Zhang PetscObjectContainer container; 1052499ec78SHong Zhang 1062499ec78SHong Zhang PetscFunctionBegin; 107899cda47SBarry Smith if (A->cmap.rstart != B->rmap.rstart || A->cmap.rend != B->rmap.rend){ 108899cda47SBarry 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); 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 */ 116899cda47SBarry 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() */ 125899cda47SBarry Smith ierr = MatMerge(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 */ 1282499ec78SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 1292499ec78SHong Zhang ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr); 130f33d1a9aSHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 1319649163dSHong Zhang ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 1324a7f9380SHong Zhang mult->MatDestroy = (*C)->ops->destroy; 133abf897f8SHong Zhang mult->MatDuplicate = (*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; 1482499ec78SHong Zhang PetscObjectContainer container; 1492499ec78SHong Zhang 1502499ec78SHong Zhang PetscFunctionBegin; 151f33d1a9aSHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 1522499ec78SHong Zhang if (container) { 1532499ec78SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 15416bf3c38SBarry Smith } else { 15516bf3c38SBarry Smith SETERRQ(PETSC_ERR_PLIB,"Container does not exit"); 1562499ec78SHong Zhang } 1572499ec78SHong Zhang 1582499ec78SHong Zhang seq = &mult->B_seq; 1592499ec78SHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1602499ec78SHong Zhang mult->B_seq = *seq; 1612499ec78SHong Zhang 1622499ec78SHong Zhang seq = &mult->A_loc; 1632499ec78SHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1642499ec78SHong Zhang mult->A_loc = *seq; 1652499ec78SHong Zhang 1662499ec78SHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 1672499ec78SHong Zhang 1682499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 169899cda47SBarry Smith ierr = MatMerge(A->comm,mult->C_seq,B->cmap.n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 1702499ec78SHong Zhang PetscFunctionReturn(0); 1712499ec78SHong Zhang } 1729123193aSHong Zhang 1739123193aSHong Zhang #undef __FUNCT__ 1749123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 1759123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1769123193aSHong Zhang { 1779123193aSHong Zhang PetscErrorCode ierr; 1789123193aSHong Zhang 1799123193aSHong Zhang PetscFunctionBegin; 1809123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 1819123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 1829123193aSHong Zhang } 1839123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 1849123193aSHong Zhang PetscFunctionReturn(0); 1859123193aSHong Zhang } 1869123193aSHong Zhang 1874324174eSBarry Smith typedef struct { 1884324174eSBarry Smith Mat workB; 1894324174eSBarry Smith PetscScalar *rvalues,*svalues; 1904324174eSBarry Smith MPI_Request *rwaits,*swaits; 1914324174eSBarry Smith } MPIAIJ_MPIDense; 1924324174eSBarry Smith 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; 211f32d5d43SBarry Smith PetscInt nz = aij->B->cmap.n; 2124324174eSBarry Smith PetscObjectContainer 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; 2179123193aSHong Zhang 2189123193aSHong Zhang PetscFunctionBegin; 2199123193aSHong Zhang ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,0.0,C); 220f32d5d43SBarry Smith 2214324174eSBarry Smith ierr = PetscObjectContainerCreate(A->comm,&cont);CHKERRQ(ierr); 2224324174eSBarry Smith ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 2234324174eSBarry Smith ierr = PetscObjectContainerSetPointer(cont,contents);CHKERRQ(ierr); 2244324174eSBarry Smith ierr = PetscObjectContainerSetUserDestroy(cont,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 2254324174eSBarry Smith 226f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 2274324174eSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap.N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 2284324174eSBarry Smith 2294324174eSBarry Smith /* Create work arrays needed */ 2304324174eSBarry Smith ierr = PetscMalloc4(B->cmap.N*from->starts[from->n],PetscScalar,&contents->rvalues, 2314324174eSBarry Smith B->cmap.N*to->starts[to->n],PetscScalar,&contents->svalues, 2324324174eSBarry Smith from->n,MPI_Request,&contents->rwaits, 2334324174eSBarry Smith to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 2344324174eSBarry Smith 2354324174eSBarry Smith ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)cont);CHKERRQ(ierr); 2364324174eSBarry Smith ierr = PetscObjectContainerDestroy(cont);CHKERRQ(ierr); 2379123193aSHong Zhang PetscFunctionReturn(0); 2389123193aSHong Zhang } 2399123193aSHong Zhang 240f32d5d43SBarry Smith /* 241*2636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 242*2636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 243f32d5d43SBarry Smith */ 2444324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 245f32d5d43SBarry Smith { 246f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 247f32d5d43SBarry Smith PetscErrorCode ierr; 248f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 249f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 250f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 251f32d5d43SBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 252f32d5d43SBarry Smith PetscInt i,j,k; 253f32d5d43SBarry Smith PetscMPIInt *sindices,*sstarts,*sprocs,*rindices,*rstarts,*rprocs,nrecvs; 254f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 255f32d5d43SBarry Smith MPI_Comm comm = A->comm; 256f32d5d43SBarry Smith PetscMPIInt tag = ctx->tag,ncols = B->cmap.N, nrows = aij->B->cmap.n,imdex,nrowsB = B->rmap.n; 257f32d5d43SBarry Smith MPI_Status status; 2584324174eSBarry Smith MPIAIJ_MPIDense *contents; 2594324174eSBarry Smith PetscObjectContainer cont; 2604324174eSBarry Smith Mat workB; 261f32d5d43SBarry Smith 262f32d5d43SBarry Smith PetscFunctionBegin; 2634324174eSBarry Smith ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&cont);CHKERRQ(ierr); 2644324174eSBarry Smith ierr = PetscObjectContainerGetPointer(cont,(void**)&contents);CHKERRQ(ierr); 2654324174eSBarry Smith 2664324174eSBarry Smith workB = *outworkB = contents->workB; 267f32d5d43SBarry 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); 268f32d5d43SBarry Smith sindices = to->indices; 269f32d5d43SBarry Smith sstarts = to->starts; 270f32d5d43SBarry Smith sprocs = to->procs; 2714324174eSBarry Smith swaits = contents->swaits; 2724324174eSBarry Smith svalues = contents->svalues; 273f32d5d43SBarry Smith 274f32d5d43SBarry Smith rindices = from->indices; 275f32d5d43SBarry Smith rstarts = from->starts; 276f32d5d43SBarry Smith rprocs = from->procs; 2774324174eSBarry Smith rwaits = contents->rwaits; 2784324174eSBarry Smith rvalues = contents->rvalues; 279f32d5d43SBarry Smith 280f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 281f32d5d43SBarry Smith ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 282f32d5d43SBarry Smith 283f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 284f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 285f32d5d43SBarry Smith } 286f32d5d43SBarry Smith 287f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 288f32d5d43SBarry Smith /* pack a message at a time */ 289f32d5d43SBarry Smith CHKMEMQ; 290f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 291f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 292f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 293*2636ff26SBarry Smith } 294*2636ff26SBarry Smith } 295f32d5d43SBarry Smith CHKMEMQ; 296f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 297f32d5d43SBarry Smith } 298*2636ff26SBarry Smith 299f32d5d43SBarry Smith nrecvs = from->n; 300f32d5d43SBarry Smith while (nrecvs) { 301f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 302f32d5d43SBarry Smith nrecvs--; 303f32d5d43SBarry Smith /* unpack a message at a time */ 304f32d5d43SBarry Smith CHKMEMQ; 305f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 306f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 307f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 308*2636ff26SBarry Smith } 309*2636ff26SBarry Smith } 310f32d5d43SBarry Smith CHKMEMQ; 311f32d5d43SBarry Smith } 312f32d5d43SBarry Smith if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr)} 313f32d5d43SBarry Smith 314f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 315f32d5d43SBarry Smith ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 316f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 317f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 318f32d5d43SBarry Smith PetscFunctionReturn(0); 319f32d5d43SBarry Smith } 320f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 321f32d5d43SBarry Smith 3229123193aSHong Zhang #undef __FUNCT__ 3239123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 3249123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 3259123193aSHong Zhang { 3269123193aSHong Zhang PetscErrorCode ierr; 327f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 328f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 329f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 330f32d5d43SBarry Smith Mat workB; 3319123193aSHong Zhang 3329123193aSHong Zhang PetscFunctionBegin; 3339123193aSHong Zhang 334f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 335f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 336f32d5d43SBarry Smith 337f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 3384324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 339f32d5d43SBarry Smith 340f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 341f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 3429123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3439123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3449123193aSHong Zhang PetscFunctionReturn(0); 3459123193aSHong Zhang } 3469123193aSHong Zhang 347