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 #include "../src/mat/impls/dense/mpi/mpidense.h" 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 15 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 16 { 17 PetscErrorCode ierr; 18 19 PetscFunctionBegin; 20 if (scall == MAT_INITIAL_MATRIX){ 21 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */ 22 } else if (scall == MAT_REUSE_MATRIX){ 23 ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr); 24 } else { 25 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 26 } 27 PetscFunctionReturn(0); 28 } 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI" 32 PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr) 33 { 34 PetscErrorCode ierr; 35 Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)ptr; 36 37 PetscFunctionBegin; 38 ierr = PetscFree2(mult->startsj,mult->startsj_r);CHKERRQ(ierr); 39 ierr = PetscFree(mult->bufa);CHKERRQ(ierr); 40 if (mult->isrowa){ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);} 41 if (mult->isrowb){ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);} 42 if (mult->iscolb){ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);} 43 if (mult->C_seq){ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);} 44 if (mult->A_loc){ierr = MatDestroy(mult->A_loc);CHKERRQ(ierr); } 45 if (mult->B_seq){ierr = MatDestroy(mult->B_seq);CHKERRQ(ierr);} 46 if (mult->B_loc){ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr);} 47 if (mult->B_oth){ierr = MatDestroy(mult->B_oth);CHKERRQ(ierr);} 48 ierr = PetscFree(mult->abi);CHKERRQ(ierr); 49 ierr = PetscFree(mult->abj);CHKERRQ(ierr); 50 ierr = PetscFree(mult);CHKERRQ(ierr); 51 PetscFunctionReturn(0); 52 } 53 54 extern PetscErrorCode MatDestroy_AIJ(Mat); 55 #undef __FUNCT__ 56 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 57 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 58 { 59 PetscErrorCode ierr; 60 PetscContainer container; 61 Mat_MatMatMultMPI *mult=PETSC_NULL; 62 63 PetscFunctionBegin; 64 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 65 if (container) { 66 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 67 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 68 A->ops->destroy = mult->MatDestroy; 69 ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 70 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 71 ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 72 PetscFunctionReturn(0); 73 } 74 75 #undef __FUNCT__ 76 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 77 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 78 { 79 PetscErrorCode ierr; 80 Mat_MatMatMultMPI *mult; 81 PetscContainer container; 82 83 PetscFunctionBegin; 84 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 85 if (container) { 86 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 87 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 88 /* Note: the container is not duplicated, because it requires deep copying of 89 several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()). 90 These data sets are only used for repeated calling of MatMatMultNumeric(). 91 *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */ 92 ierr = (*mult->MatDuplicate)(A,op,M);CHKERRQ(ierr); 93 (*M)->ops->destroy = mult->MatDestroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 94 (*M)->ops->duplicate = mult->MatDuplicate; /* = MatDuplicate_MPIAIJ */ 95 PetscFunctionReturn(0); 96 } 97 98 #undef __FUNCT__ 99 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 100 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 101 { 102 PetscErrorCode ierr; 103 PetscInt start,end; 104 Mat_MatMatMultMPI *mult; 105 PetscContainer container; 106 107 PetscFunctionBegin; 108 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 109 SETERRQ4(PETSC_COMM_SELF,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); 110 } 111 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 112 113 /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 114 ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr); 115 116 /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 117 start = A->rmap->rstart; end = A->rmap->rend; 118 ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr); 119 ierr = MatGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr); 120 121 /* compute C_seq = A_seq * B_seq */ 122 ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr); 123 124 /* create mpi matrix C by concatinating C_seq */ 125 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 126 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 127 128 /* attach the supporting struct to C for reuse of symbolic C */ 129 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 130 ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 131 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 132 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 133 mult->MatDestroy = (*C)->ops->destroy; 134 mult->MatDuplicate = (*C)->ops->duplicate; 135 136 (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 137 (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 138 PetscFunctionReturn(0); 139 } 140 141 /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 142 #undef __FUNCT__ 143 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 144 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 145 { 146 PetscErrorCode ierr; 147 Mat *seq; 148 Mat_MatMatMultMPI *mult; 149 PetscContainer container; 150 151 PetscFunctionBegin; 152 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 153 if (container) { 154 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 155 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 156 157 seq = &mult->B_seq; 158 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 159 mult->B_seq = *seq; 160 161 seq = &mult->A_loc; 162 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 163 mult->A_loc = *seq; 164 165 ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 166 167 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 168 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 169 PetscFunctionReturn(0); 170 } 171 172 #undef __FUNCT__ 173 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 174 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 175 { 176 PetscErrorCode ierr; 177 178 PetscFunctionBegin; 179 if (scall == MAT_INITIAL_MATRIX){ 180 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 181 } 182 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 183 PetscFunctionReturn(0); 184 } 185 186 typedef struct { 187 Mat workB; 188 PetscScalar *rvalues,*svalues; 189 MPI_Request *rwaits,*swaits; 190 } MPIAIJ_MPIDense; 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 194 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 195 { 196 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 197 PetscErrorCode ierr; 198 199 PetscFunctionBegin; 200 if (contents->workB) {ierr = MatDestroy(contents->workB);CHKERRQ(ierr);} 201 ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 202 ierr = PetscFree(contents);CHKERRQ(ierr); 203 PetscFunctionReturn(0); 204 } 205 206 #undef __FUNCT__ 207 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 208 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 209 { 210 PetscErrorCode ierr; 211 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 212 PetscInt nz = aij->B->cmap->n; 213 PetscContainer cont; 214 MPIAIJ_MPIDense *contents; 215 VecScatter ctx = aij->Mvctx; 216 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 217 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 218 PetscInt m=A->rmap->n,n=B->cmap->n; 219 220 PetscFunctionBegin; 221 ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 222 ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 223 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 224 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 225 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 226 227 ierr = PetscContainerCreate(((PetscObject)A)->comm,&cont);CHKERRQ(ierr); 228 ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 229 ierr = PetscContainerSetPointer(cont,contents);CHKERRQ(ierr); 230 ierr = PetscContainerSetUserDestroy(cont,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 231 232 /* Create work matrix used to store off processor rows of B needed for local product */ 233 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 234 235 /* Create work arrays needed */ 236 ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 237 B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 238 from->n,MPI_Request,&contents->rwaits, 239 to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 240 241 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)cont);CHKERRQ(ierr); 242 ierr = PetscContainerDestroy(cont);CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "MatMPIDenseScatter" 248 /* 249 Performs an efficient scatter on the rows of B needed by this process; this is 250 a modification of the VecScatterBegin_() routines. 251 */ 252 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 253 { 254 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 255 PetscErrorCode ierr; 256 PetscScalar *b,*w,*svalues,*rvalues; 257 VecScatter ctx = aij->Mvctx; 258 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 259 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 260 PetscInt i,j,k; 261 PetscInt *sindices,*sstarts,*rindices,*rstarts; 262 PetscMPIInt *sprocs,*rprocs,nrecvs; 263 MPI_Request *swaits,*rwaits; 264 MPI_Comm comm = ((PetscObject)A)->comm; 265 PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 266 MPI_Status status; 267 MPIAIJ_MPIDense *contents; 268 PetscContainer cont; 269 Mat workB; 270 271 PetscFunctionBegin; 272 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&cont);CHKERRQ(ierr); 273 ierr = PetscContainerGetPointer(cont,(void**)&contents);CHKERRQ(ierr); 274 275 workB = *outworkB = contents->workB; 276 if (nrows != workB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n); 277 sindices = to->indices; 278 sstarts = to->starts; 279 sprocs = to->procs; 280 swaits = contents->swaits; 281 svalues = contents->svalues; 282 283 rindices = from->indices; 284 rstarts = from->starts; 285 rprocs = from->procs; 286 rwaits = contents->rwaits; 287 rvalues = contents->rvalues; 288 289 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 290 ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 291 292 for (i=0; i<from->n; i++) { 293 ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 294 } 295 296 for (i=0; i<to->n; i++) { 297 /* pack a message at a time */ 298 CHKMEMQ; 299 for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 300 for (k=0; k<ncols; k++) { 301 svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 302 } 303 } 304 CHKMEMQ; 305 ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 306 } 307 308 nrecvs = from->n; 309 while (nrecvs) { 310 ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 311 nrecvs--; 312 /* unpack a message at a time */ 313 CHKMEMQ; 314 for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 315 for (k=0; k<ncols; k++) { 316 w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 317 } 318 } 319 CHKMEMQ; 320 } 321 if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 322 323 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 324 ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 325 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 326 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 327 PetscFunctionReturn(0); 328 } 329 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 330 331 #undef __FUNCT__ 332 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 333 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 334 { 335 PetscErrorCode ierr; 336 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 337 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 338 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 339 Mat workB; 340 341 PetscFunctionBegin; 342 343 /* diagonal block of A times all local rows of B*/ 344 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 345 346 /* get off processor parts of B needed to complete the product */ 347 ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 348 349 /* off-diagonal block of A times nonlocal rows of B */ 350 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 351 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 353 PetscFunctionReturn(0); 354 } 355 356