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