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:   PetscFree(mult->startsj);
 38:   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:   PetscFree(mult->abi);
 48:   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_PLIB,"Container does not exit");
 68:   }
 69:   A->ops->destroy = mult->MatDestroy;
 70:   PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);
 71:   (*A->ops->destroy)(A);
 72:   PetscObjectContainerDestroy(container);
 73:   return(0);
 74: }

 78: PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) {
 79:   PetscErrorCode       ierr;
 80:   Mat_MatMatMultMPI    *mult;
 81:   PetscObjectContainer container;

 84:   PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);
 85:   if (container) {
 86:     PetscObjectContainerGetPointer(container,(void **)&mult);
 87:   } else {
 88:     SETERRQ(PETSC_ERR_PLIB,"Container does not exit");
 89:   }
 90:   (*mult->MatDuplicate)(A,op,M);
 91:   (*M)->ops->destroy   = mult->MatDestroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
 92:   (*M)->ops->duplicate = mult->MatDuplicate; /* =MatDuplicate_ MPIAIJ */
 93:   return(0);
 94: }

 98: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
 99: {
100:   PetscErrorCode       ierr;
101:   PetscInt             start,end;
102:   Mat_MatMatMultMPI    *mult;
103:   PetscObjectContainer container;
104: 
106:   if (A->cmap.rstart != B->rmap.rstart || A->cmap.rend != B->rmap.rend){
107:     SETERRQ4(PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap.rstart,A->cmap.rend,B->rmap.rstart,B->rmap.rend);
108:   }
109:   PetscNew(Mat_MatMatMultMPI,&mult);

111:   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
112:   MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);

114:   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
115:   start = A->rmap.rstart; end = A->rmap.rend;
116:   ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);
117:   MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);

119:   /* compute C_seq = A_seq * B_seq */
120:   MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);

122:   /* create mpi matrix C by concatinating C_seq */
123:   PetscObjectReference((PetscObject)mult->C_seq); /* prevent C_seq being destroyed by MatMerge() */
124:   MatMerge(A->comm,mult->C_seq,B->cmap.n,MAT_INITIAL_MATRIX,C);
125: 
126:   /* attach the supporting struct to C for reuse of symbolic C */
127:   PetscObjectContainerCreate(PETSC_COMM_SELF,&container);
128:   PetscObjectContainerSetPointer(container,mult);
129:   PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);
130:   PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);
131:   mult->MatDestroy   = (*C)->ops->destroy;
132:   mult->MatDuplicate = (*C)->ops->duplicate;

134:   (*C)->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
135:   (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
136:   return(0);
137: }

139: /* This routine is called ONLY in the case of reusing previously computed symbolic C */
142: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
143: {
144:   PetscErrorCode       ierr;
145:   Mat                  *seq;
146:   Mat_MatMatMultMPI    *mult;
147:   PetscObjectContainer container;

150:   PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);
151:   if (container) {
152:     PetscObjectContainerGetPointer(container,(void **)&mult);
153:   } else {
154:     SETERRQ(PETSC_ERR_PLIB,"Container does not exit");
155:   }

157:   seq = &mult->B_seq;
158:   MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);
159:   mult->B_seq = *seq;
160: 
161:   seq = &mult->A_loc;
162:   MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);
163:   mult->A_loc = *seq;

165:   MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);

167:   PetscObjectReference((PetscObject)mult->C_seq);
168:   MatMerge(A->comm,mult->C_seq,B->cmap.n,MAT_REUSE_MATRIX,&C);
169:   return(0);
170: }