xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision f33d1a9a8e1cc388a7592d5011fcd885988d41e4)
12499ec78SHong Zhang /*
22499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
32499ec78SHong Zhang           C = A * B
42499ec78SHong Zhang */
52499ec78SHong Zhang 
62499ec78SHong Zhang #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
72499ec78SHong Zhang #include "src/mat/utils/freespace.h"
82499ec78SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
92499ec78SHong Zhang #include "petscbt.h"
102499ec78SHong Zhang 
112499ec78SHong Zhang #undef __FUNCT__
122499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
132499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
142499ec78SHong Zhang {
152499ec78SHong Zhang   PetscErrorCode ierr;
162499ec78SHong Zhang 
172499ec78SHong Zhang   PetscFunctionBegin;
182499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
192499ec78SHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
202499ec78SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
212499ec78SHong Zhang     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
222499ec78SHong Zhang   } else {
232499ec78SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
242499ec78SHong Zhang   }
252499ec78SHong Zhang   PetscFunctionReturn(0);
262499ec78SHong Zhang }
272499ec78SHong Zhang 
28*f33d1a9aSHong Zhang #undef __FUNCT__
29*f33d1a9aSHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_Mat_MatMatMultMPI"
30*f33d1a9aSHong Zhang PetscErrorCode PetscObjectContainerDestroy_Mat_MatMatMultMPI(PetscObject obj)
31*f33d1a9aSHong Zhang {
32*f33d1a9aSHong Zhang   PetscErrorCode       ierr;
33*f33d1a9aSHong Zhang   Mat                  A=(Mat)obj;
34*f33d1a9aSHong Zhang   PetscObjectContainer container;
35*f33d1a9aSHong Zhang   Mat_MatMatMultMPI    *mult=PETSC_NULL;
36*f33d1a9aSHong Zhang 
37*f33d1a9aSHong Zhang   PetscFunctionBegin;
38*f33d1a9aSHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
39*f33d1a9aSHong Zhang   if (container) {
40*f33d1a9aSHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
41*f33d1a9aSHong Zhang     if (mult->startsj){ierr = PetscFree(mult->startsj);CHKERRQ(ierr);}
42*f33d1a9aSHong Zhang     if (mult->bufa){ierr = PetscFree(mult->bufa);CHKERRQ(ierr);}
43*f33d1a9aSHong Zhang     if (mult->isrowa){ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);}
44*f33d1a9aSHong Zhang     if (mult->isrowb){ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);}
45*f33d1a9aSHong Zhang     if (mult->iscolb){ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);}
46*f33d1a9aSHong Zhang     if (mult->C_seq){ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);}
47*f33d1a9aSHong Zhang     if (mult->A_loc){ierr = MatDestroy(mult->A_loc);CHKERRQ(ierr); }
48*f33d1a9aSHong Zhang     if (mult->B_seq){ierr = MatDestroy(mult->B_seq);CHKERRQ(ierr);}
49*f33d1a9aSHong Zhang     if (mult->B_loc){ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr);}
50*f33d1a9aSHong Zhang     if (mult->B_oth){ierr = MatDestroy(mult->B_oth);CHKERRQ(ierr);}
51*f33d1a9aSHong Zhang     if (mult->abi){ierr = PetscFree(mult->abi);CHKERRQ(ierr);}
52*f33d1a9aSHong Zhang     if (mult->abj){ierr = PetscFree(mult->abj);CHKERRQ(ierr);}
53*f33d1a9aSHong Zhang     ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
54*f33d1a9aSHong Zhang     ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
55*f33d1a9aSHong Zhang     ierr = PetscFree(mult);CHKERRQ(ierr);
56*f33d1a9aSHong Zhang   }
57*f33d1a9aSHong Zhang   PetscFunctionReturn(0);
58*f33d1a9aSHong Zhang }
59*f33d1a9aSHong Zhang 
602499ec78SHong Zhang EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
612499ec78SHong Zhang #undef __FUNCT__
622499ec78SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
632499ec78SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
642499ec78SHong Zhang {
652499ec78SHong Zhang   PetscErrorCode       ierr;
662499ec78SHong Zhang 
672499ec78SHong Zhang   PetscFunctionBegin;
68*f33d1a9aSHong Zhang   ierr = PetscObjectContainerDestroy_Mat_MatMatMultMPI((PetscObject)A);CHKERRQ(ierr);
692499ec78SHong Zhang   ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
702499ec78SHong Zhang   PetscFunctionReturn(0);
712499ec78SHong Zhang }
722499ec78SHong Zhang 
732499ec78SHong Zhang #undef __FUNCT__
742499ec78SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
752499ec78SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
762499ec78SHong Zhang {
772499ec78SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*b=(Mat_MPIAIJ*)B->data;
782499ec78SHong Zhang   PetscErrorCode       ierr;
792499ec78SHong Zhang   PetscInt             start,end;
802499ec78SHong Zhang   Mat_MatMatMultMPI    *mult;
812499ec78SHong Zhang   PetscObjectContainer container;
822499ec78SHong Zhang 
832499ec78SHong Zhang   PetscFunctionBegin;
842499ec78SHong Zhang   if (a->cstart != b->rstart || a->cend != b->rend){
852499ec78SHong Zhang     SETERRQ4(PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",a->cstart,a->cend,b->rstart,b->rend);
862499ec78SHong Zhang   }
872499ec78SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
88*f33d1a9aSHong Zhang   ierr = PetscMemzero(mult,sizeof(Mat_MatMatMultMPI));CHKERRQ(ierr);
892499ec78SHong Zhang 
902499ec78SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
912499ec78SHong Zhang   ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr);
922499ec78SHong Zhang 
932499ec78SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
942499ec78SHong Zhang   start = a->rstart; end = a->rend;
952499ec78SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
962499ec78SHong Zhang   ierr = MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr);
972499ec78SHong Zhang 
982499ec78SHong Zhang   /* compute C_seq = A_seq * B_seq */
992499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
1002499ec78SHong Zhang 
1012499ec78SHong Zhang   /* create mpi matrix C by concatinating C_seq */
1022499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
1032499ec78SHong Zhang   ierr = MatMerge(A->comm,mult->C_seq,B->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
1042499ec78SHong Zhang 
1052499ec78SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
1062499ec78SHong Zhang   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
1072499ec78SHong Zhang   ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr);
108*f33d1a9aSHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
1092499ec78SHong Zhang 
1102499ec78SHong Zhang   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
1112499ec78SHong Zhang 
1122499ec78SHong Zhang   PetscFunctionReturn(0);
1132499ec78SHong Zhang }
1142499ec78SHong Zhang 
1152499ec78SHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
1162499ec78SHong Zhang #undef __FUNCT__
1172499ec78SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
1182499ec78SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
1192499ec78SHong Zhang {
1202499ec78SHong Zhang   PetscErrorCode       ierr;
1212499ec78SHong Zhang   Mat                  *seq;
1222499ec78SHong Zhang   Mat_MatMatMultMPI    *mult;
1232499ec78SHong Zhang   PetscObjectContainer container;
1242499ec78SHong Zhang 
1252499ec78SHong Zhang   PetscFunctionBegin;
126*f33d1a9aSHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
1272499ec78SHong Zhang   if (container) {
1282499ec78SHong Zhang     ierr  = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
1292499ec78SHong Zhang   }
1302499ec78SHong Zhang 
1312499ec78SHong Zhang   seq = &mult->B_seq;
1322499ec78SHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
1332499ec78SHong Zhang   mult->B_seq = *seq;
1342499ec78SHong Zhang 
1352499ec78SHong Zhang   seq = &mult->A_loc;
1362499ec78SHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
1372499ec78SHong Zhang   mult->A_loc = *seq;
1382499ec78SHong Zhang 
1392499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
1402499ec78SHong Zhang 
1412499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
1422499ec78SHong Zhang   ierr = MatMerge(A->comm,mult->C_seq,B->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
1432499ec78SHong Zhang 
1442499ec78SHong Zhang   PetscFunctionReturn(0);
1452499ec78SHong Zhang }
146