xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision d44834fb5aa6ab55eaa2e18f09dfd00e40aa4639)
1 /*
2   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
3           C = A * B
4 */
5 
6 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
7 #include "src/mat/utils/freespace.h"
8 #include "src/mat/impls/aij/mpi/mpiaij.h"
9 #include "petscbt.h"
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
13 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
14 {
15   PetscErrorCode ierr;
16 
17   PetscFunctionBegin;
18   if (scall == MAT_INITIAL_MATRIX){
19     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
20   } else if (scall == MAT_REUSE_MATRIX){
21     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
22   } else {
23     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
24   }
25   PetscFunctionReturn(0);
26 }
27 
28 #undef __FUNCT__
29 #define __FUNCT__ "PetscObjectContainerDestroy_Mat_MatMatMultMPI"
30 PetscErrorCode PetscObjectContainerDestroy_Mat_MatMatMultMPI(PetscObject obj)
31 {
32   PetscErrorCode       ierr;
33   Mat                  A=(Mat)obj;
34   PetscObjectContainer container;
35   Mat_MatMatMultMPI    *mult=PETSC_NULL;
36 
37   PetscFunctionBegin;
38   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
39   if (container) {
40     ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
41     if (mult->startsj){ierr = PetscFree(mult->startsj);CHKERRQ(ierr);}
42     if (mult->bufa){ierr = PetscFree(mult->bufa);CHKERRQ(ierr);}
43     if (mult->isrowa){ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);}
44     if (mult->isrowb){ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);}
45     if (mult->iscolb){ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);}
46     if (mult->C_seq){ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);}
47     if (mult->A_loc){ierr = MatDestroy(mult->A_loc);CHKERRQ(ierr); }
48     if (mult->B_seq){ierr = MatDestroy(mult->B_seq);CHKERRQ(ierr);}
49     if (mult->B_loc){ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr);}
50     if (mult->B_oth){ierr = MatDestroy(mult->B_oth);CHKERRQ(ierr);}
51     if (mult->abi){ierr = PetscFree(mult->abi);CHKERRQ(ierr);}
52     if (mult->abj){ierr = PetscFree(mult->abj);CHKERRQ(ierr);}
53     ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
54     ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
55     ierr = PetscFree(mult);CHKERRQ(ierr);
56   }
57   PetscFunctionReturn(0);
58 }
59 
60 EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
61 #undef __FUNCT__
62 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
63 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
64 {
65   PetscErrorCode       ierr;
66 
67   PetscFunctionBegin;
68   ierr = PetscObjectContainerDestroy_Mat_MatMatMultMPI((PetscObject)A);CHKERRQ(ierr);
69   ierr = MatDestroy_MPIAIJ(A);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   ierr = PetscMemzero(mult,sizeof(Mat_MatMatMultMPI));CHKERRQ(ierr);
89 
90   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
91   ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr);
92 
93   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
94   start = a->rstart; end = a->rend;
95   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
96   ierr = MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr);
97 
98   /* compute C_seq = A_seq * B_seq */
99   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
100 
101   /* create mpi matrix C by concatinating C_seq */
102   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
103   ierr = MatMerge(A->comm,mult->C_seq,B->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
104 
105   /* attach the supporting struct to C for reuse of symbolic C */
106   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
107   ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr);
108   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
109 
110   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
111 
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 
144   PetscFunctionReturn(0);
145 }
146