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