Actual source code: mpimatmatmult.c

  1: /*
  2:   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
  3:           C = A * B
  4: */

 6:  #include src/mat/impls/aij/seq/aij.h
 7:  #include src/mat/utils/freespace.h
 8:  #include src/mat/impls/aij/mpi/mpiaij.h
 9:  #include petscbt.h

 13: PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
 14: {

 18:   if (scall == MAT_INITIAL_MATRIX){
 19:     MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);/* numeric product is computed as well */
 20:   } else if (scall == MAT_REUSE_MATRIX){
 21:     MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);
 22:   } else {
 23:     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
 24:   }
 25:   return(0);
 26: }

 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;

 38:   PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);
 39:   if (container) {
 40:     PetscObjectContainerGetPointer(container,(void **)&mult);
 41:     if (mult->startsj){PetscFree(mult->startsj);}
 42:     if (mult->bufa){PetscFree(mult->bufa);}
 43:     if (mult->isrowa){ISDestroy(mult->isrowa);}
 44:     if (mult->isrowb){ISDestroy(mult->isrowb);}
 45:     if (mult->iscolb){ISDestroy(mult->iscolb);}
 46:     if (mult->C_seq){MatDestroy(mult->C_seq);}
 47:     if (mult->A_loc){MatDestroy(mult->A_loc); }
 48:     if (mult->B_seq){MatDestroy(mult->B_seq);}
 49:     if (mult->B_loc){MatDestroy(mult->B_loc);}
 50:     if (mult->B_oth){MatDestroy(mult->B_oth);}
 51:     if (mult->abi){PetscFree(mult->abi);}
 52:     if (mult->abj){PetscFree(mult->abj);}
 53:     PetscObjectContainerDestroy(container);
 54:     PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);
 55:     PetscFree(mult);
 56:   }
 57:   return(0);
 58: }

 60: EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
 63: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
 64: {
 65:   PetscErrorCode       ierr;

 68:   PetscObjectContainerDestroy_Mat_MatMatMultMPI((PetscObject)A);
 69:   MatDestroy_MPIAIJ(A);
 70:   return(0);
 71: }

 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: 
 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:   PetscNew(Mat_MatMatMultMPI,&mult);
 88:   PetscMemzero(mult,sizeof(Mat_MatMatMultMPI));

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

 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:   ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);
 96:   MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);

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

101:   /* create mpi matrix C by concatinating C_seq */
102:   PetscObjectReference((PetscObject)mult->C_seq); /* prevent C_seq being destroyed by MatMerge() */
103:   MatMerge(A->comm,mult->C_seq,B->n,MAT_INITIAL_MATRIX,C);
104: 
105:   /* attach the supporting struct to C for reuse of symbolic C */
106:   PetscObjectContainerCreate(PETSC_COMM_SELF,&container);
107:   PetscObjectContainerSetPointer(container,mult);
108:   PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);

110:   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
111: 
112:   return(0);
113: }

115: /* This routine is called ONLY in the case of reusing previously computed symbolic C */
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;

126:   PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);
127:   if (container) {
128:     PetscObjectContainerGetPointer(container,(void **)&mult);
129:   }

131:   seq = &mult->B_seq;
132:   MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);
133:   mult->B_seq = *seq;
134: 
135:   seq = &mult->A_loc;
136:   MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);
137:   mult->A_loc = *seq;

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

141:   PetscObjectReference((PetscObject)mult->C_seq);
142:   MatMerge(A->comm,mult->C_seq,B->n,MAT_REUSE_MATRIX,&C);

144:   return(0);
145: }