1 2 /* 3 Defines matrix-matrix product routines for pairs of MPIAIJ matrices 4 C = A * B 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 #include <../src/mat/impls/dense/mpi/mpidense.h> 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 14 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 15 { 16 PetscErrorCode ierr; 17 18 PetscFunctionBegin; 19 if (scall == MAT_INITIAL_MATRIX){ 20 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 21 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well in old version! */ 22 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 23 } 24 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 25 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 26 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 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 = ISDestroy(&mult->isrowa);CHKERRQ(ierr); 39 ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr); 40 ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr); 41 ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr); 42 ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr); 43 ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr); 44 ierr = PetscFree(mult);CHKERRQ(ierr); 45 PetscFunctionReturn(0); 46 } 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_new" 50 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_new(Mat A) 51 { 52 PetscErrorCode ierr; 53 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 54 Mat_PtAPMPI *ptap=a->ptap; 55 56 PetscFunctionBegin; 57 ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr); 58 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 59 ierr = MatDestroy(&ptap->B_loc);CHKERRQ(ierr); 60 ierr = MatDestroy(&ptap->B_oth);CHKERRQ(ierr); 61 ierr = PetscFree(ptap->abi);CHKERRQ(ierr); 62 ierr = PetscFree(ptap->abj);CHKERRQ(ierr); 63 ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 64 ierr = ptap->destroy(A);CHKERRQ(ierr); 65 ierr = PetscFree(ptap);CHKERRQ(ierr); 66 PetscFunctionReturn(0); 67 } 68 69 #undef __FUNCT__ 70 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 71 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 72 { 73 PetscErrorCode ierr; 74 PetscContainer container; 75 Mat_MatMatMultMPI *mult=PETSC_NULL; 76 77 PetscFunctionBegin; 78 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 79 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 80 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 81 A->ops->destroy = mult->destroy; 82 A->ops->duplicate = mult->duplicate; 83 if (A->ops->destroy) { 84 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 85 } 86 ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 87 PetscFunctionReturn(0); 88 } 89 90 #undef __FUNCT__ 91 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 92 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 93 { 94 PetscErrorCode ierr; 95 Mat_MatMatMultMPI *mult; 96 PetscContainer container; 97 98 PetscFunctionBegin; 99 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 100 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 101 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 102 /* Note: the container is not duplicated, because it requires deep copying of 103 several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()). 104 These data sets are only used for repeated calling of MatMatMultNumeric(). 105 *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */ 106 ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr); 107 (*M)->ops->destroy = mult->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 108 (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */ 109 PetscFunctionReturn(0); 110 } 111 112 #undef __FUNCT__ 113 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_new" 114 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_new(Mat A, MatDuplicateOption op, Mat *M) 115 { 116 PetscErrorCode ierr; 117 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 118 Mat_PtAPMPI *ptap=a->ptap; 119 120 PetscFunctionBegin; 121 ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 122 (*M)->ops->destroy = ptap->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */ 123 (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */ 124 PetscFunctionReturn(0); 125 } 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_new" 129 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_new(Mat A,Mat P,Mat C) 130 { 131 PetscErrorCode ierr; 132 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 133 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 134 Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 135 PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 136 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 137 Mat_SeqAIJ *p_loc,*p_oth; 138 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 139 PetscScalar *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca; 140 PetscInt cm=C->rmap->n,anz,pnz; 141 Mat_PtAPMPI *ptap=c->ptap; 142 PetscInt *api,*apj,*apJ,cnz,i,j,k,row; 143 PetscInt rstart=C->rmap->rstart,cstart=C->cmap->rstart; 144 PetscInt cdnz,conz,k0,k1; 145 146 PetscFunctionBegin; 147 /* 1) get P_oth = ptap->B_oth and P_loc = ptap->B_loc */ 148 /*-----------------------------------------------------*/ 149 if (ptap->reuse == MAT_INITIAL_MATRIX){ 150 ptap->reuse = MAT_REUSE_MATRIX; 151 } else { /* update numerical values of P_oth and P_loc */ 152 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->B_oth);CHKERRQ(ierr); 153 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->B_loc);CHKERRQ(ierr); 154 } 155 156 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 157 /*----------------------------------------------------------*/ 158 /* get data from symbolic products */ 159 p_loc = (Mat_SeqAIJ*)(ptap->B_loc)->data; 160 p_oth = (Mat_SeqAIJ*)(ptap->B_oth)->data; 161 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 162 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 163 164 /* get apa for storing dense row A[i,:]*P */ 165 apa = ptap->apa; 166 167 for (i=0; i<cm; i++) { 168 /* diagonal portion of A */ 169 anz = adi[i+1] - adi[i]; 170 adj = ad->j + adi[i]; 171 ada = ad->a + adi[i]; 172 for (j=0; j<anz; j++) { 173 row = adj[j]; 174 pnz = pi_loc[row+1] - pi_loc[row]; 175 pj = pj_loc + pi_loc[row]; 176 pa = pa_loc + pi_loc[row]; 177 178 /* perform dense axpy */ 179 valtmp = ada[j]; 180 for (k=0; k<pnz; k++){ 181 apa[pj[k]] += valtmp*pa[k]; 182 } 183 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 184 } 185 186 /* off-diagonal portion of A */ 187 anz = aoi[i+1] - aoi[i]; 188 aoj = ao->j + aoi[i]; 189 aoa = ao->a + aoi[i]; 190 for (j=0; j<anz; j++) { 191 row = aoj[j]; 192 pnz = pi_oth[row+1] - pi_oth[row]; 193 pj = pj_oth + pi_oth[row]; 194 pa = pa_oth + pi_oth[row]; 195 196 /* perform dense axpy */ 197 valtmp = aoa[j]; 198 for (k=0; k<pnz; k++){ 199 apa[pj[k]] += valtmp*pa[k]; 200 } 201 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 202 } 203 204 /* set values in C */ 205 row = rstart + i; 206 api = ptap->abi; 207 apj = ptap->abj; 208 apJ = apj + api[i]; 209 cnz = api[i+1] - api[i]; 210 cdnz = cd->i[i+1] - cd->i[i]; 211 conz = co->i[i+1] - co->i[i]; 212 213 /* 1st off-diagoanl part of C */ 214 ca = coa + co->i[i]; 215 k = 0; 216 for (k0=0; k0<conz; k0++){ 217 if (apJ[k] >= cstart) break; 218 ca[k0] = apa[apJ[k]]; 219 apa[apJ[k]] = 0.0; 220 k++; 221 } 222 223 /* diagonal part of C */ 224 ca = cda + cd->i[i]; 225 for (k1=0; k1<cdnz; k1++){ 226 ca[k1] = apa[apJ[k]]; 227 apa[apJ[k]] = 0.0; 228 k++; 229 } 230 231 /* 2nd off-diagoanl part of C */ 232 ca = coa + co->i[i]; 233 for (; k0<conz; k0++){ 234 ca[k0] = apa[apJ[k]]; 235 apa[apJ[k]] = 0.0; 236 k++; 237 } 238 } 239 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 240 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 241 PetscFunctionReturn(0); 242 } 243 244 #undef __FUNCT__ 245 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_new" 246 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_new(Mat A,Mat P,PetscReal fill,Mat *C) 247 { 248 PetscErrorCode ierr; 249 MPI_Comm comm=((PetscObject)A)->comm; 250 PetscMPIInt size,rank; 251 Mat B_mpi; 252 Mat_PtAPMPI *ptap; 253 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 254 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 255 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 256 Mat_SeqAIJ *p_loc,*p_oth; 257 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 258 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 259 PetscInt nlnk,*lnk,i,pnz,row; 260 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 261 PetscBT lnkbt; 262 PetscInt nzi; 263 PetscInt *api,*apj,*Jptr,apnz,nspacedouble=0,j; 264 PetscScalar *apa; 265 266 PetscFunctionBegin; 267 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 268 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 269 270 /* create struct Mat_PtAPMPI and attached it to C later */ 271 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 272 ptap->abnz_max = 0; 273 ptap->reuse = MAT_INITIAL_MATRIX; 274 275 /* malloc apa to store dense row A[i,:]*P */ 276 ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 277 ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 278 ptap->apa = apa; 279 280 281 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 282 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->B_oth);CHKERRQ(ierr); 283 /* get P_loc by taking all local rows of P */ 284 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->B_loc);CHKERRQ(ierr); 285 286 p_loc = (Mat_SeqAIJ*)(ptap->B_loc)->data; 287 p_oth = (Mat_SeqAIJ*)(ptap->B_oth)->data; 288 pi_loc = p_loc->i; pj_loc = p_loc->j; 289 pi_oth = p_oth->i; pj_oth = p_oth->j; 290 291 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 292 /*-------------------------------------------------------------------*/ 293 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 294 ptap->abi = api; 295 api[0] = 0; 296 297 /* create and initialize a linked list */ 298 nlnk = pN+1; 299 ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 300 301 /* Initial FreeSpace size is fill*nnz(A) */ 302 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr); 303 current_space = free_space; 304 305 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 306 for (i=0;i<am;i++) { 307 apnz = 0; 308 /* diagonal portion of A */ 309 nzi = adi[i+1] - adi[i]; 310 for (j=0; j<nzi; j++){ 311 row = *adj++; 312 pnz = pi_loc[row+1] - pi_loc[row]; 313 Jptr = pj_loc + pi_loc[row]; 314 /* add non-zero cols of P into the sorted linked list lnk */ 315 ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 316 apnz += nlnk; 317 } 318 /* off-diagonal portion of A */ 319 nzi = aoi[i+1] - aoi[i]; 320 for (j=0; j<nzi; j++){ 321 row = *aoj++; 322 pnz = pi_oth[row+1] - pi_oth[row]; 323 Jptr = pj_oth + pi_oth[row]; 324 ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 325 apnz += nlnk; 326 } 327 328 api[i+1] = api[i] + apnz; 329 if (ptap->abnz_max < apnz) ptap->abnz_max = apnz; 330 331 /* if free space is not available, double the total space in the list */ 332 if (current_space->local_remaining<apnz) { 333 ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 334 nspacedouble++; 335 } 336 337 /* Copy data into free space, then initialize lnk */ 338 ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 339 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 340 current_space->array += apnz; 341 current_space->local_used += apnz; 342 current_space->local_remaining -= apnz; 343 } 344 345 /* Allocate space for apj, initialize apj, and */ 346 /* destroy list of free space and other temporary array(s) */ 347 ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->abj);CHKERRQ(ierr); 348 apj = ptap->abj; 349 ierr = PetscFreeSpaceContiguous(&free_space,ptap->abj);CHKERRQ(ierr); 350 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 351 352 /* create and assemble symbolic parallel matrix B_mpi */ 353 /*----------------------------------------------------*/ 354 ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr); 355 ierr = MatSetSizes(B_mpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 356 ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr); 357 ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr); 358 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 359 ierr = MatSetBlockSize(B_mpi,1);CHKERRQ(ierr); 360 for (i=0; i<am; i++){ 361 row = i + rstart; 362 apnz = api[i+1] - api[i]; 363 ierr = MatSetValues(B_mpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 364 apj += apnz; 365 } 366 ierr = MatAssemblyBegin(B_mpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 367 ierr = MatAssemblyEnd(B_mpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 368 369 ptap->destroy = B_mpi->ops->destroy; 370 ptap->duplicate = B_mpi->ops->duplicate; 371 B_mpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_new; 372 B_mpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_new; 373 B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_new; 374 375 /* attach the supporting struct to B_mpi for reuse */ 376 c = (Mat_MPIAIJ*)B_mpi->data; 377 c->ptap = ptap; 378 379 *C = B_mpi; 380 PetscFunctionReturn(0); 381 } 382 383 #undef __FUNCT__ 384 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 385 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 386 { 387 PetscErrorCode ierr; 388 Mat_MatMatMultMPI *mult; 389 PetscContainer container; 390 Mat AB,*seq; 391 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 392 PetscInt *idx,i,start,ncols,nzA,nzB,*cmap,imark; 393 PetscBool matmatmult_new=PETSC_FALSE; 394 395 PetscFunctionBegin; 396 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 397 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); 398 } 399 ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_new",&matmatmult_new,PETSC_NULL);CHKERRQ(ierr); 400 if (matmatmult_new){ 401 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_new(A,B,fill,C);;CHKERRQ(ierr); 402 PetscFunctionReturn(0); 403 } 404 405 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 406 407 /* get isrowb: nonzero col of A */ 408 start = A->cmap->rstart; 409 cmap = a->garray; 410 nzA = a->A->cmap->n; 411 nzB = a->B->cmap->n; 412 ierr = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr); 413 ncols = 0; 414 for (i=0; i<nzB; i++) { /* row < local row index */ 415 if (cmap[i] < start) idx[ncols++] = cmap[i]; 416 else break; 417 } 418 imark = i; 419 for (i=0; i<nzA; i++) idx[ncols++] = start + i; /* local rows */ 420 for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */ 421 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr); 422 ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr); 423 424 /* get isrowa: all local rows of A */ 425 ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr); 426 427 /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */ 428 /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 429 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 430 mult->B_seq = *seq; 431 ierr = PetscFree(seq);CHKERRQ(ierr); 432 433 /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 434 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 435 mult->A_loc = *seq; 436 ierr = PetscFree(seq);CHKERRQ(ierr); 437 438 /* compute C_seq = A_seq * B_seq */ 439 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr); 440 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 441 442 /* create mpi matrix C by concatinating C_seq */ 443 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 444 ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr); 445 ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr); 446 447 /* attach the supporting struct to C for reuse of symbolic C */ 448 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 449 ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 450 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 451 ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 452 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 453 mult->destroy = AB->ops->destroy; 454 mult->duplicate = AB->ops->duplicate; 455 AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 456 AB->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 457 AB->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 458 AB->ops->matmult = MatMatMult_MPIAIJ_MPIAIJ; 459 *C = AB; 460 PetscFunctionReturn(0); 461 } 462 463 /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 464 #undef __FUNCT__ 465 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 466 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 467 { 468 PetscErrorCode ierr; 469 Mat *seq; 470 Mat_MatMatMultMPI *mult; 471 PetscContainer container; 472 473 PetscFunctionBegin; 474 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 475 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 476 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 477 seq = &mult->B_seq; 478 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 479 mult->B_seq = *seq; 480 481 seq = &mult->A_loc; 482 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 483 mult->A_loc = *seq; 484 485 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 486 487 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 488 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 492 #undef __FUNCT__ 493 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 494 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 495 { 496 PetscErrorCode ierr; 497 498 PetscFunctionBegin; 499 if (scall == MAT_INITIAL_MATRIX){ 500 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 501 } 502 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 503 PetscFunctionReturn(0); 504 } 505 506 typedef struct { 507 Mat workB; 508 PetscScalar *rvalues,*svalues; 509 MPI_Request *rwaits,*swaits; 510 } MPIAIJ_MPIDense; 511 512 #undef __FUNCT__ 513 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 514 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 515 { 516 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 517 PetscErrorCode ierr; 518 519 PetscFunctionBegin; 520 ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 521 ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 522 ierr = PetscFree(contents);CHKERRQ(ierr); 523 PetscFunctionReturn(0); 524 } 525 526 #undef __FUNCT__ 527 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 528 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 529 { 530 PetscErrorCode ierr; 531 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 532 PetscInt nz = aij->B->cmap->n; 533 PetscContainer container; 534 MPIAIJ_MPIDense *contents; 535 VecScatter ctx = aij->Mvctx; 536 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 537 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 538 PetscInt m=A->rmap->n,n=B->cmap->n; 539 540 PetscFunctionBegin; 541 ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 542 ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 543 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 544 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 545 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 546 (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 547 548 ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 549 /* Create work matrix used to store off processor rows of B needed for local product */ 550 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 551 /* Create work arrays needed */ 552 ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 553 B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 554 from->n,MPI_Request,&contents->rwaits, 555 to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 556 557 ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 558 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 559 ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 560 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 561 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 562 PetscFunctionReturn(0); 563 } 564 565 #undef __FUNCT__ 566 #define __FUNCT__ "MatMPIDenseScatter" 567 /* 568 Performs an efficient scatter on the rows of B needed by this process; this is 569 a modification of the VecScatterBegin_() routines. 570 */ 571 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 572 { 573 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 574 PetscErrorCode ierr; 575 PetscScalar *b,*w,*svalues,*rvalues; 576 VecScatter ctx = aij->Mvctx; 577 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 578 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 579 PetscInt i,j,k; 580 PetscInt *sindices,*sstarts,*rindices,*rstarts; 581 PetscMPIInt *sprocs,*rprocs,nrecvs; 582 MPI_Request *swaits,*rwaits; 583 MPI_Comm comm = ((PetscObject)A)->comm; 584 PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 585 MPI_Status status; 586 MPIAIJ_MPIDense *contents; 587 PetscContainer container; 588 Mat workB; 589 590 PetscFunctionBegin; 591 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 592 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 593 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 594 595 workB = *outworkB = contents->workB; 596 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); 597 sindices = to->indices; 598 sstarts = to->starts; 599 sprocs = to->procs; 600 swaits = contents->swaits; 601 svalues = contents->svalues; 602 603 rindices = from->indices; 604 rstarts = from->starts; 605 rprocs = from->procs; 606 rwaits = contents->rwaits; 607 rvalues = contents->rvalues; 608 609 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 610 ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 611 612 for (i=0; i<from->n; i++) { 613 ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 614 } 615 616 for (i=0; i<to->n; i++) { 617 /* pack a message at a time */ 618 CHKMEMQ; 619 for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 620 for (k=0; k<ncols; k++) { 621 svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 622 } 623 } 624 CHKMEMQ; 625 ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 626 } 627 628 nrecvs = from->n; 629 while (nrecvs) { 630 ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 631 nrecvs--; 632 /* unpack a message at a time */ 633 CHKMEMQ; 634 for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 635 for (k=0; k<ncols; k++) { 636 w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 637 } 638 } 639 CHKMEMQ; 640 } 641 if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 642 643 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 644 ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 645 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 646 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 647 PetscFunctionReturn(0); 648 } 649 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 650 651 #undef __FUNCT__ 652 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 653 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 654 { 655 PetscErrorCode ierr; 656 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 657 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 658 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 659 Mat workB; 660 661 PetscFunctionBegin; 662 663 /* diagonal block of A times all local rows of B*/ 664 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 665 666 /* get off processor parts of B needed to complete the product */ 667 ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 668 669 /* off-diagonal block of A times nonlocal rows of B */ 670 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 671 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 672 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 673 PetscFunctionReturn(0); 674 } 675 676