xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 2499ec7893278572238d73a6b308a0d02f10e2c6)
1*2499ec78SHong Zhang /*
2*2499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
3*2499ec78SHong Zhang           C = A * B
4*2499ec78SHong Zhang */
5*2499ec78SHong Zhang 
6*2499ec78SHong Zhang #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
7*2499ec78SHong Zhang #include "src/mat/utils/freespace.h"
8*2499ec78SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
9*2499ec78SHong Zhang #include "petscbt.h"
10*2499ec78SHong Zhang 
11*2499ec78SHong Zhang #undef __FUNCT__
12*2499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
13*2499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
14*2499ec78SHong Zhang {
15*2499ec78SHong Zhang   PetscErrorCode ierr;
16*2499ec78SHong Zhang 
17*2499ec78SHong Zhang   PetscFunctionBegin;
18*2499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
19*2499ec78SHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
20*2499ec78SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
21*2499ec78SHong Zhang     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
22*2499ec78SHong Zhang   } else {
23*2499ec78SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
24*2499ec78SHong Zhang   }
25*2499ec78SHong Zhang   PetscFunctionReturn(0);
26*2499ec78SHong Zhang }
27*2499ec78SHong Zhang 
28*2499ec78SHong Zhang EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
29*2499ec78SHong Zhang #undef __FUNCT__
30*2499ec78SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
31*2499ec78SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
32*2499ec78SHong Zhang {
33*2499ec78SHong Zhang   PetscErrorCode       ierr;
34*2499ec78SHong Zhang   Mat_MatMatMultMPI    *mult;
35*2499ec78SHong Zhang   PetscObjectContainer container;
36*2499ec78SHong Zhang 
37*2499ec78SHong Zhang   PetscFunctionBegin;
38*2499ec78SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
39*2499ec78SHong Zhang   if (container) {
40*2499ec78SHong Zhang     ierr  = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
41*2499ec78SHong Zhang     ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);
42*2499ec78SHong Zhang     ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);
43*2499ec78SHong Zhang     ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);
44*2499ec78SHong Zhang     ierr = MatDestroy(mult->A_loc);CHKERRQ(ierr);
45*2499ec78SHong Zhang     ierr = MatDestroy(mult->B_seq);CHKERRQ(ierr);
46*2499ec78SHong Zhang     ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);
47*2499ec78SHong Zhang 
48*2499ec78SHong Zhang     ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
49*2499ec78SHong Zhang     ierr = PetscObjectCompose((PetscObject)A,"MatMatMultMPI",0);CHKERRQ(ierr);
50*2499ec78SHong Zhang   }
51*2499ec78SHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
52*2499ec78SHong Zhang   ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
53*2499ec78SHong Zhang 
54*2499ec78SHong Zhang   PetscFunctionReturn(0);
55*2499ec78SHong Zhang }
56*2499ec78SHong Zhang 
57*2499ec78SHong Zhang #undef __FUNCT__
58*2499ec78SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
59*2499ec78SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
60*2499ec78SHong Zhang {
61*2499ec78SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*b=(Mat_MPIAIJ*)B->data;
62*2499ec78SHong Zhang   PetscErrorCode       ierr;
63*2499ec78SHong Zhang   PetscInt             start,end;
64*2499ec78SHong Zhang   Mat_MatMatMultMPI    *mult;
65*2499ec78SHong Zhang   PetscObjectContainer container;
66*2499ec78SHong Zhang 
67*2499ec78SHong Zhang   PetscFunctionBegin;
68*2499ec78SHong Zhang   if (a->cstart != b->rstart || a->cend != b->rend){
69*2499ec78SHong Zhang     SETERRQ4(PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",a->cstart,a->cend,b->rstart,b->rend);
70*2499ec78SHong Zhang   }
71*2499ec78SHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
72*2499ec78SHong Zhang 
73*2499ec78SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
74*2499ec78SHong Zhang   ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr);
75*2499ec78SHong Zhang 
76*2499ec78SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
77*2499ec78SHong Zhang   start = a->rstart; end = a->rend;
78*2499ec78SHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
79*2499ec78SHong Zhang   ierr = MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr);
80*2499ec78SHong Zhang 
81*2499ec78SHong Zhang   /* compute C_seq = A_seq * B_seq */
82*2499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
83*2499ec78SHong Zhang 
84*2499ec78SHong Zhang   /* create mpi matrix C by concatinating C_seq */
85*2499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
86*2499ec78SHong Zhang   ierr = MatMerge(A->comm,mult->C_seq,B->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
87*2499ec78SHong Zhang 
88*2499ec78SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
89*2499ec78SHong Zhang   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
90*2499ec78SHong Zhang   ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr);
91*2499ec78SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
92*2499ec78SHong Zhang 
93*2499ec78SHong Zhang   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
94*2499ec78SHong Zhang 
95*2499ec78SHong Zhang   PetscFunctionReturn(0);
96*2499ec78SHong Zhang }
97*2499ec78SHong Zhang 
98*2499ec78SHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
99*2499ec78SHong Zhang #undef __FUNCT__
100*2499ec78SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
101*2499ec78SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
102*2499ec78SHong Zhang {
103*2499ec78SHong Zhang   PetscErrorCode       ierr;
104*2499ec78SHong Zhang   Mat                  *seq;
105*2499ec78SHong Zhang   Mat_MatMatMultMPI    *mult;
106*2499ec78SHong Zhang   PetscObjectContainer container;
107*2499ec78SHong Zhang 
108*2499ec78SHong Zhang   PetscFunctionBegin;
109*2499ec78SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
110*2499ec78SHong Zhang   if (container) {
111*2499ec78SHong Zhang     ierr  = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
112*2499ec78SHong Zhang   }
113*2499ec78SHong Zhang 
114*2499ec78SHong Zhang   seq = &mult->B_seq;
115*2499ec78SHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
116*2499ec78SHong Zhang   mult->B_seq = *seq;
117*2499ec78SHong Zhang 
118*2499ec78SHong Zhang   seq = &mult->A_loc;
119*2499ec78SHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
120*2499ec78SHong Zhang   mult->A_loc = *seq;
121*2499ec78SHong Zhang 
122*2499ec78SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
123*2499ec78SHong Zhang 
124*2499ec78SHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
125*2499ec78SHong Zhang   ierr = MatMerge(A->comm,mult->C_seq,B->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
126*2499ec78SHong Zhang 
127*2499ec78SHong Zhang   PetscFunctionReturn(0);
128*2499ec78SHong Zhang }
129