xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 4ae313f42201736c22765e1d4ea08232bbac4af6)
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"
11*4ae313f4SHong 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 
1879123193aSHong Zhang #undef __FUNCT__
1889123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
1899123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
1909123193aSHong Zhang {
1919123193aSHong Zhang   PetscErrorCode ierr;
1929123193aSHong Zhang 
1939123193aSHong Zhang   PetscFunctionBegin;
1949123193aSHong Zhang   ierr = MatMatMultSymbolic_MPIDense_MPIDense(A,B,0.0,C);
1959123193aSHong Zhang   PetscFunctionReturn(0);
1969123193aSHong Zhang }
1979123193aSHong Zhang 
1989123193aSHong Zhang #undef __FUNCT__
1999123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
2009123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
2019123193aSHong Zhang {
2029123193aSHong Zhang   PetscErrorCode ierr;
2039123193aSHong Zhang   PetscScalar    *b,*c,*barray,*carray;
2049123193aSHong Zhang   PetscInt       cm=C->rmap.n, bN=B->cmap.N, bm=B->rmap.n, col;
2059123193aSHong Zhang   Vec            vb,vc;
2069123193aSHong Zhang 
2079123193aSHong Zhang   PetscFunctionBegin;
2089123193aSHong Zhang   if (!C->rmap.N || !bN) PetscFunctionReturn(0);
2099123193aSHong Zhang   ierr = VecCreateMPI(B->comm,bm,PETSC_DECIDE,&vb);CHKERRQ(ierr);
2109123193aSHong Zhang   ierr = VecCreateMPI(A->comm,cm,PETSC_DECIDE,&vc);CHKERRQ(ierr);
2119123193aSHong Zhang 
2129123193aSHong Zhang   ierr = MatGetArray(B,&barray);CHKERRQ(ierr);
2139123193aSHong Zhang   ierr = MatGetArray(C,&carray);CHKERRQ(ierr);
2149123193aSHong Zhang   b = barray; c = carray;
2159123193aSHong Zhang   for (col=0; col<bN; col++){
2169123193aSHong Zhang     ierr = VecPlaceArray(vb,b);CHKERRQ(ierr);
2179123193aSHong Zhang     ierr = VecPlaceArray(vc,c);CHKERRQ(ierr);
2189123193aSHong Zhang     ierr = MatMult(A,vb,vc);CHKERRQ(ierr);
2199123193aSHong Zhang     b += bm; c += cm;
2209123193aSHong Zhang     ierr = VecResetArray(vb);CHKERRQ(ierr);
2219123193aSHong Zhang     ierr = VecResetArray(vc);CHKERRQ(ierr);
2229123193aSHong Zhang   }
2239123193aSHong Zhang   ierr = MatRestoreArray(B,&barray);CHKERRQ(ierr);
2249123193aSHong Zhang   ierr = MatRestoreArray(C,&carray);CHKERRQ(ierr);
2259123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2269123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2279123193aSHong Zhang   ierr = VecDestroy(vb);CHKERRQ(ierr);
2289123193aSHong Zhang   ierr = VecDestroy(vc);CHKERRQ(ierr);
2299123193aSHong Zhang   PetscFunctionReturn(0);
2309123193aSHong Zhang }
2319123193aSHong Zhang 
232