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); 64bf0cc555SLisandro 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; 67bf0cc555SLisandro Dalcin A->ops->duplicate = mult->duplicate; 68bf0cc555SLisandro Dalcin if (A->ops->destroy) { 69992edd73SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 70bf0cc555SLisandro Dalcin } 71bf0cc555SLisandro 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); 85bf0cc555SLisandro 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 Mat_MatMatMultMPI *mult; 103776b82aeSLisandro Dalcin PetscContainer container; 104d5821860SHong Zhang Mat AB,*seq; 105d5821860SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 106d5821860SHong Zhang PetscInt *idx,i,start,ncols,nzA,nzB,*cmap,imark; 1072499ec78SHong Zhang 1082499ec78SHong Zhang PetscFunctionBegin; 109d0f46423SBarry Smith if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 110e32f2f54SBarry 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); 1112499ec78SHong Zhang } 1122499ec78SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 1132499ec78SHong Zhang 114d5821860SHong Zhang /* get isrowb: nonzero col of A */ 115d5821860SHong Zhang start = A->cmap->rstart; 116d5821860SHong Zhang cmap = a->garray; 117d5821860SHong Zhang nzA = a->A->cmap->n; 118d5821860SHong Zhang nzB = a->B->cmap->n; 119d5821860SHong Zhang ierr = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr); 120d5821860SHong Zhang ncols = 0; 121d5821860SHong Zhang for (i=0; i<nzB; i++) { /* row < local row index */ 122d5821860SHong Zhang if (cmap[i] < start) idx[ncols++] = cmap[i]; 123d5821860SHong Zhang else break; 124d5821860SHong Zhang } 125d5821860SHong Zhang imark = i; 126d5821860SHong Zhang for (i=0; i<nzA; i++) idx[ncols++] = start + i; /* local rows */ 127d5821860SHong Zhang for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */ 128d5821860SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr); 129d5821860SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr); 130d5821860SHong Zhang 131d5821860SHong Zhang /* get isrowa: all local rows of A */ 132d5821860SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr); 133d5821860SHong Zhang 134d5821860SHong Zhang /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */ 1352499ec78SHong Zhang /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 136d5821860SHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 137d5821860SHong Zhang mult->B_seq = *seq; 138d5821860SHong Zhang ierr = PetscFree(seq);CHKERRQ(ierr); 1392499ec78SHong Zhang 1402499ec78SHong Zhang /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 141d5821860SHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 142d5821860SHong Zhang mult->A_loc = *seq; 143d5821860SHong Zhang ierr = PetscFree(seq);CHKERRQ(ierr); 1442499ec78SHong Zhang 1452499ec78SHong Zhang /* compute C_seq = A_seq * B_seq */ 146*9b8102ccSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr); 147*9b8102ccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 1482499ec78SHong Zhang 1492499ec78SHong Zhang /* create mpi matrix C by concatinating C_seq */ 1502499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 151*9b8102ccSHong Zhang ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr); 152*9b8102ccSHong Zhang ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr); 1532499ec78SHong Zhang 1542499ec78SHong Zhang /* attach the supporting struct to C for reuse of symbolic C */ 155776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 156776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 157776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 158d5821860SHong Zhang ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 159b160f555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 160d5821860SHong Zhang mult->destroy = AB->ops->destroy; 161d5821860SHong Zhang mult->duplicate = AB->ops->duplicate; 162d5821860SHong Zhang AB->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 163d5821860SHong Zhang AB->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 1648cdbd757SHong Zhang AB->ops->matmult = MatMatMult_MPIAIJ_MPIAIJ; 165d5821860SHong Zhang *C = AB; 1662499ec78SHong Zhang PetscFunctionReturn(0); 1672499ec78SHong Zhang } 1682499ec78SHong Zhang 1692499ec78SHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 1702499ec78SHong Zhang #undef __FUNCT__ 1712499ec78SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 1722499ec78SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 1732499ec78SHong Zhang { 1742499ec78SHong Zhang PetscErrorCode ierr; 1752499ec78SHong Zhang Mat *seq; 1762499ec78SHong Zhang Mat_MatMatMultMPI *mult; 177776b82aeSLisandro Dalcin PetscContainer container; 1782499ec78SHong Zhang 1792499ec78SHong Zhang PetscFunctionBegin; 180f33d1a9aSHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 181bf0cc555SLisandro Dalcin if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 182776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 1832499ec78SHong Zhang seq = &mult->B_seq; 1842499ec78SHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1852499ec78SHong Zhang mult->B_seq = *seq; 1862499ec78SHong Zhang 1872499ec78SHong Zhang seq = &mult->A_loc; 1882499ec78SHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1892499ec78SHong Zhang mult->A_loc = *seq; 1902499ec78SHong Zhang 1912499ec78SHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 1922499ec78SHong Zhang 1932499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 194d0f46423SBarry Smith ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 1952499ec78SHong Zhang PetscFunctionReturn(0); 1962499ec78SHong Zhang } 1979123193aSHong Zhang 1989123193aSHong Zhang #undef __FUNCT__ 1999123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 2009123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 2019123193aSHong Zhang { 2029123193aSHong Zhang PetscErrorCode ierr; 2039123193aSHong Zhang 2049123193aSHong Zhang PetscFunctionBegin; 2059123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2069123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 2079123193aSHong Zhang } 2089123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 2099123193aSHong Zhang PetscFunctionReturn(0); 2109123193aSHong Zhang } 2119123193aSHong Zhang 2124324174eSBarry Smith typedef struct { 2134324174eSBarry Smith Mat workB; 2144324174eSBarry Smith PetscScalar *rvalues,*svalues; 2154324174eSBarry Smith MPI_Request *rwaits,*swaits; 2164324174eSBarry Smith } MPIAIJ_MPIDense; 2174324174eSBarry Smith 2187af9e4fdSHong Zhang #undef __FUNCT__ 2197af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 2204324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 2214324174eSBarry Smith { 2224324174eSBarry Smith MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 2234324174eSBarry Smith PetscErrorCode ierr; 2244324174eSBarry Smith 2254324174eSBarry Smith PetscFunctionBegin; 2266bf464f9SBarry Smith ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 2274324174eSBarry Smith ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 2284324174eSBarry Smith ierr = PetscFree(contents);CHKERRQ(ierr); 2294324174eSBarry Smith PetscFunctionReturn(0); 2304324174eSBarry Smith } 2314324174eSBarry Smith 2329123193aSHong Zhang #undef __FUNCT__ 2339123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 2349123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 2359123193aSHong Zhang { 2369123193aSHong Zhang PetscErrorCode ierr; 237f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 238d0f46423SBarry Smith PetscInt nz = aij->B->cmap->n; 239bf0cc555SLisandro Dalcin PetscContainer container; 2404324174eSBarry Smith MPIAIJ_MPIDense *contents; 2414324174eSBarry Smith VecScatter ctx = aij->Mvctx; 2424324174eSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 2434324174eSBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 244d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 2459123193aSHong Zhang 2469123193aSHong Zhang PetscFunctionBegin; 247cb2480eaSBarry Smith ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 248d0f46423SBarry Smith ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 249cb2480eaSBarry Smith ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 250cb2480eaSBarry Smith ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 251cb2480eaSBarry Smith ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2528cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 253f32d5d43SBarry Smith 2544324174eSBarry Smith ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 255f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 256d0f46423SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 2574324174eSBarry Smith /* Create work arrays needed */ 258d0f46423SBarry Smith ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 259d0f46423SBarry Smith B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 2604324174eSBarry Smith from->n,MPI_Request,&contents->rwaits, 2614324174eSBarry Smith to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 2624324174eSBarry Smith 263bf0cc555SLisandro Dalcin ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 264bf0cc555SLisandro Dalcin ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 265bf0cc555SLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 266bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 267bf0cc555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 2689123193aSHong Zhang PetscFunctionReturn(0); 2699123193aSHong Zhang } 2709123193aSHong Zhang 2717af9e4fdSHong Zhang #undef __FUNCT__ 2727af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter" 273f32d5d43SBarry Smith /* 2742636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 2752636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 276f32d5d43SBarry Smith */ 2774324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 278f32d5d43SBarry Smith { 279f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 280f32d5d43SBarry Smith PetscErrorCode ierr; 281f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 282f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 283f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 284f32d5d43SBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 285f32d5d43SBarry Smith PetscInt i,j,k; 286aa5bb8c0SSatish Balay PetscInt *sindices,*sstarts,*rindices,*rstarts; 287aa5bb8c0SSatish Balay PetscMPIInt *sprocs,*rprocs,nrecvs; 288f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 2897adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)A)->comm; 290d0f46423SBarry Smith PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 291f32d5d43SBarry Smith MPI_Status status; 2924324174eSBarry Smith MPIAIJ_MPIDense *contents; 293bf0cc555SLisandro Dalcin PetscContainer container; 2944324174eSBarry Smith Mat workB; 295f32d5d43SBarry Smith 296f32d5d43SBarry Smith PetscFunctionBegin; 297bf0cc555SLisandro Dalcin ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 298bf0cc555SLisandro Dalcin if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 299bf0cc555SLisandro Dalcin ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 3004324174eSBarry Smith 3014324174eSBarry Smith workB = *outworkB = contents->workB; 302e32f2f54SBarry 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); 303f32d5d43SBarry Smith sindices = to->indices; 304f32d5d43SBarry Smith sstarts = to->starts; 305f32d5d43SBarry Smith sprocs = to->procs; 3064324174eSBarry Smith swaits = contents->swaits; 3074324174eSBarry Smith svalues = contents->svalues; 308f32d5d43SBarry Smith 309f32d5d43SBarry Smith rindices = from->indices; 310f32d5d43SBarry Smith rstarts = from->starts; 311f32d5d43SBarry Smith rprocs = from->procs; 3124324174eSBarry Smith rwaits = contents->rwaits; 3134324174eSBarry Smith rvalues = contents->rvalues; 314f32d5d43SBarry Smith 315f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 316f32d5d43SBarry Smith ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 317f32d5d43SBarry Smith 318f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 319f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 320f32d5d43SBarry Smith } 321f32d5d43SBarry Smith 322f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 323f32d5d43SBarry Smith /* pack a message at a time */ 324f32d5d43SBarry Smith CHKMEMQ; 325f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 326f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 327f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 3282636ff26SBarry Smith } 3292636ff26SBarry Smith } 330f32d5d43SBarry Smith CHKMEMQ; 331f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 332f32d5d43SBarry Smith } 3332636ff26SBarry Smith 334f32d5d43SBarry Smith nrecvs = from->n; 335f32d5d43SBarry Smith while (nrecvs) { 336f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 337f32d5d43SBarry Smith nrecvs--; 338f32d5d43SBarry Smith /* unpack a message at a time */ 339f32d5d43SBarry Smith CHKMEMQ; 340f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 341f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 342f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 3432636ff26SBarry Smith } 3442636ff26SBarry Smith } 345f32d5d43SBarry Smith CHKMEMQ; 346f32d5d43SBarry Smith } 347cb9801acSJed Brown if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 348f32d5d43SBarry Smith 349f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 350f32d5d43SBarry Smith ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 351f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 353f32d5d43SBarry Smith PetscFunctionReturn(0); 354f32d5d43SBarry Smith } 355f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 356f32d5d43SBarry Smith 3579123193aSHong Zhang #undef __FUNCT__ 3589123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 3599123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 3609123193aSHong Zhang { 3619123193aSHong Zhang PetscErrorCode ierr; 362f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 363f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 364f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 365f32d5d43SBarry Smith Mat workB; 3669123193aSHong Zhang 3679123193aSHong Zhang PetscFunctionBegin; 3689123193aSHong Zhang 369f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 370f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 371f32d5d43SBarry Smith 372f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 3734324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 374f32d5d43SBarry Smith 375f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 376f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 3779123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3789123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3799123193aSHong Zhang PetscFunctionReturn(0); 3809123193aSHong Zhang } 3819123193aSHong Zhang 382