1*be1d678aSKris Buschelman #define PETSCMAT_DLL 2*be1d678aSKris 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" 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 { 242499ec78SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 252499ec78SHong Zhang } 262499ec78SHong Zhang PetscFunctionReturn(0); 272499ec78SHong Zhang } 282499ec78SHong Zhang 29f33d1a9aSHong Zhang #undef __FUNCT__ 30f33d1a9aSHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_Mat_MatMatMultMPI" 319649163dSHong Zhang PetscErrorCode PetscObjectContainerDestroy_Mat_MatMatMultMPI(void *ptr) 32f33d1a9aSHong Zhang { 33f33d1a9aSHong Zhang PetscErrorCode ierr; 349649163dSHong Zhang Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)ptr; 35f33d1a9aSHong Zhang 36f33d1a9aSHong Zhang PetscFunctionBegin; 37f33d1a9aSHong Zhang if (mult->startsj){ierr = PetscFree(mult->startsj);CHKERRQ(ierr);} 38f33d1a9aSHong Zhang if (mult->bufa){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);} 47f33d1a9aSHong Zhang if (mult->abi){ierr = PetscFree(mult->abi);CHKERRQ(ierr);} 48f33d1a9aSHong Zhang if (mult->abj){ierr = PetscFree(mult->abj);CHKERRQ(ierr);} 49f33d1a9aSHong Zhang ierr = PetscFree(mult);CHKERRQ(ierr); 50f33d1a9aSHong Zhang PetscFunctionReturn(0); 51f33d1a9aSHong Zhang } 52f33d1a9aSHong Zhang 534a7f9380SHong Zhang 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; 599649163dSHong Zhang PetscObjectContainer 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) { 659649163dSHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 669649163dSHong Zhang } else { 679649163dSHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 689649163dSHong Zhang } 699649163dSHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 706c6e00e0SHong Zhang ierr = (*mult->MatDestroy)(A);CHKERRQ(ierr); 719649163dSHong Zhang ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 722499ec78SHong Zhang PetscFunctionReturn(0); 732499ec78SHong Zhang } 742499ec78SHong Zhang 752499ec78SHong Zhang #undef __FUNCT__ 762499ec78SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 772499ec78SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 782499ec78SHong Zhang { 792499ec78SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*b=(Mat_MPIAIJ*)B->data; 802499ec78SHong Zhang PetscErrorCode ierr; 812499ec78SHong Zhang PetscInt start,end; 822499ec78SHong Zhang Mat_MatMatMultMPI *mult; 832499ec78SHong Zhang PetscObjectContainer container; 842499ec78SHong Zhang 852499ec78SHong Zhang PetscFunctionBegin; 862499ec78SHong Zhang if (a->cstart != b->rstart || a->cend != b->rend){ 872499ec78SHong Zhang SETERRQ4(PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",a->cstart,a->cend,b->rstart,b->rend); 882499ec78SHong Zhang } 892499ec78SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 902499ec78SHong Zhang 912499ec78SHong Zhang /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 922499ec78SHong Zhang ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr); 932499ec78SHong Zhang 942499ec78SHong Zhang /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 952499ec78SHong Zhang start = a->rstart; end = a->rend; 962499ec78SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr); 972499ec78SHong Zhang ierr = MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr); 982499ec78SHong Zhang 992499ec78SHong Zhang /* compute C_seq = A_seq * B_seq */ 1002499ec78SHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr); 1012499ec78SHong Zhang 1022499ec78SHong Zhang /* create mpi matrix C by concatinating C_seq */ 1032499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 1042499ec78SHong Zhang ierr = MatMerge(A->comm,mult->C_seq,B->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 1052499ec78SHong Zhang 1062499ec78SHong Zhang /* attach the supporting struct to C for reuse of symbolic C */ 1072499ec78SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 1082499ec78SHong Zhang ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr); 109f33d1a9aSHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 1109649163dSHong Zhang ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 1114a7f9380SHong Zhang mult->MatDestroy = (*C)->ops->destroy; 1122499ec78SHong Zhang 1132499ec78SHong Zhang (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1142499ec78SHong Zhang PetscFunctionReturn(0); 1152499ec78SHong Zhang } 1162499ec78SHong Zhang 1172499ec78SHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 1182499ec78SHong Zhang #undef __FUNCT__ 1192499ec78SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 1202499ec78SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 1212499ec78SHong Zhang { 1222499ec78SHong Zhang PetscErrorCode ierr; 1232499ec78SHong Zhang Mat *seq; 1242499ec78SHong Zhang Mat_MatMatMultMPI *mult; 1252499ec78SHong Zhang PetscObjectContainer container; 1262499ec78SHong Zhang 1272499ec78SHong Zhang PetscFunctionBegin; 128f33d1a9aSHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 1292499ec78SHong Zhang if (container) { 1302499ec78SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 1312499ec78SHong Zhang } 1322499ec78SHong Zhang 1332499ec78SHong Zhang seq = &mult->B_seq; 1342499ec78SHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1352499ec78SHong Zhang mult->B_seq = *seq; 1362499ec78SHong Zhang 1372499ec78SHong Zhang seq = &mult->A_loc; 1382499ec78SHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 1392499ec78SHong Zhang mult->A_loc = *seq; 1402499ec78SHong Zhang 1412499ec78SHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 1422499ec78SHong Zhang 1432499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 1442499ec78SHong Zhang ierr = MatMerge(A->comm,mult->C_seq,B->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 1452499ec78SHong Zhang PetscFunctionReturn(0); 1462499ec78SHong Zhang } 147