Actual source code: mpimatmatmult.c
1: #define PETSCMAT_DLL
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
8: #include src/mat/utils/freespace.h
9: #include src/mat/impls/aij/mpi/mpiaij.h
10: #include petscbt.h
14: PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
15: {
19: if (scall == MAT_INITIAL_MATRIX){
20: MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);/* numeric product is computed as well */
21: } else if (scall == MAT_REUSE_MATRIX){
22: MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);
23: } else {
24: SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
25: }
26: return(0);
27: }
31: PetscErrorCode PetscObjectContainerDestroy_Mat_MatMatMultMPI(void *ptr)
32: {
33: PetscErrorCode ierr;
34: Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)ptr;
37: if (mult->startsj){PetscFree(mult->startsj);}
38: if (mult->bufa){PetscFree(mult->bufa);}
39: if (mult->isrowa){ISDestroy(mult->isrowa);}
40: if (mult->isrowb){ISDestroy(mult->isrowb);}
41: if (mult->iscolb){ISDestroy(mult->iscolb);}
42: if (mult->C_seq){MatDestroy(mult->C_seq);}
43: if (mult->A_loc){MatDestroy(mult->A_loc); }
44: if (mult->B_seq){MatDestroy(mult->B_seq);}
45: if (mult->B_loc){MatDestroy(mult->B_loc);}
46: if (mult->B_oth){MatDestroy(mult->B_oth);}
47: if (mult->abi){PetscFree(mult->abi);}
48: if (mult->abj){PetscFree(mult->abj);}
49: PetscFree(mult);
50: return(0);
51: }
53: EXTERN PetscErrorCode MatDestroy_AIJ(Mat);
56: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
57: {
58: PetscErrorCode ierr;
59: PetscObjectContainer container;
60: Mat_MatMatMultMPI *mult=PETSC_NULL;
63: PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);
64: if (container) {
65: PetscObjectContainerGetPointer(container,(void **)&mult);
66: } else {
67: SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
68: }
69: PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);
70: (*mult->MatDestroy)(A);
71: PetscObjectContainerDestroy(container);
72: return(0);
73: }
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:
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: PetscNew(Mat_MatMatMultMPI,&mult);
91: /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
92: MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);
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: ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);
97: MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);
99: /* compute C_seq = A_seq * B_seq */
100: MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);
102: /* create mpi matrix C by concatinating C_seq */
103: PetscObjectReference((PetscObject)mult->C_seq); /* prevent C_seq being destroyed by MatMerge() */
104: MatMerge(A->comm,mult->C_seq,B->n,MAT_INITIAL_MATRIX,C);
105:
106: /* attach the supporting struct to C for reuse of symbolic C */
107: PetscObjectContainerCreate(PETSC_COMM_SELF,&container);
108: PetscObjectContainerSetPointer(container,mult);
109: PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);
110: PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);
111: mult->MatDestroy = (*C)->ops->destroy;
113: (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
114: return(0);
115: }
117: /* This routine is called ONLY in the case of reusing previously computed symbolic C */
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;
128: PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);
129: if (container) {
130: PetscObjectContainerGetPointer(container,(void **)&mult);
131: }
133: seq = &mult->B_seq;
134: MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);
135: mult->B_seq = *seq;
136:
137: seq = &mult->A_loc;
138: MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);
139: mult->A_loc = *seq;
141: MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);
143: PetscObjectReference((PetscObject)mult->C_seq);
144: MatMerge(A->comm,mult->C_seq,B->n,MAT_REUSE_MATRIX,&C);
145: return(0);
146: }