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 #define DEBUG_MATMATMULT 13 */ 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 17 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 18 { 19 PetscErrorCode ierr; 20 21 PetscFunctionBegin; 22 if (scall == MAT_INITIAL_MATRIX){ 23 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 24 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 25 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 26 } 27 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 28 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 29 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 30 PetscFunctionReturn(0); 31 } 32 33 #undef __FUNCT__ 34 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI" 35 PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr) 36 { 37 PetscErrorCode ierr; 38 Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)ptr; 39 40 PetscFunctionBegin; 41 ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr); 42 ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr); 43 ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr); 44 ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr); 45 ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr); 46 ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr); 47 ierr = PetscFree(mult);CHKERRQ(ierr); 48 PetscFunctionReturn(0); 49 } 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 53 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 54 { 55 PetscErrorCode ierr; 56 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 57 Mat_PtAPMPI *ptap=a->ptap; 58 59 PetscFunctionBegin; 60 ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr); 61 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 62 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 63 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 64 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 65 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 66 ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 67 ierr = ptap->destroy(A);CHKERRQ(ierr); 68 ierr = PetscFree(ptap);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32" 74 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A) 75 { 76 PetscErrorCode ierr; 77 PetscContainer container; 78 Mat_MatMatMultMPI *mult=PETSC_NULL; 79 80 PetscFunctionBegin; 81 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 82 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 83 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 84 A->ops->destroy = mult->destroy; 85 A->ops->duplicate = mult->duplicate; 86 if (A->ops->destroy) { 87 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 88 } 89 ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32" 95 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M) 96 { 97 PetscErrorCode ierr; 98 Mat_MatMatMultMPI *mult; 99 PetscContainer container; 100 101 PetscFunctionBegin; 102 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 103 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 104 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 105 /* Note: the container is not duplicated, because it requires deep copying of 106 several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()). 107 These data sets are only used for repeated calling of MatMatMultNumeric(). 108 *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */ 109 ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr); 110 (*M)->ops->destroy = mult->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 111 (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */ 112 PetscFunctionReturn(0); 113 } 114 115 #undef __FUNCT__ 116 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 117 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 118 { 119 PetscErrorCode ierr; 120 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 121 Mat_PtAPMPI *ptap=a->ptap; 122 123 PetscFunctionBegin; 124 ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 125 (*M)->ops->destroy = ptap->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */ 126 (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */ 127 PetscFunctionReturn(0); 128 } 129 130 #undef __FUNCT__ 131 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 132 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 133 { 134 PetscErrorCode ierr; 135 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 136 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 137 Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 138 PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 139 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 140 Mat_SeqAIJ *p_loc,*p_oth; 141 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 142 PetscScalar *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca; 143 PetscInt cm=C->rmap->n,anz,pnz; 144 Mat_PtAPMPI *ptap=c->ptap; 145 PetscInt *api,*apj,*apJ,i,j,k,row; 146 PetscInt rstart=C->rmap->rstart,cstart=C->cmap->rstart; 147 PetscInt cdnz,conz,k0,k1; 148 #if defined(DEBUG_MATMATMULT) 149 PetscMPIInt rank=a->rank; 150 #endif 151 152 PetscFunctionBegin; 153 #if defined(DEBUG_MATMATMULT) 154 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank); 155 #endif 156 157 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 158 /*-----------------------------------------------------*/ 159 /* update numerical values of P_oth and P_loc */ 160 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 161 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 162 #if defined(DEBUG_MATMATMULT) 163 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 164 #endif 165 166 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 167 /*----------------------------------------------------------*/ 168 /* get data from symbolic products */ 169 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 170 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 171 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 172 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 173 174 /* get apa for storing dense row A[i,:]*P */ 175 apa = ptap->apa; 176 177 for (i=0; i<cm; i++) { 178 /* diagonal portion of A */ 179 anz = adi[i+1] - adi[i]; 180 adj = ad->j + adi[i]; 181 ada = ad->a + adi[i]; 182 for (j=0; j<anz; j++) { 183 row = adj[j]; 184 pnz = pi_loc[row+1] - pi_loc[row]; 185 pj = pj_loc + pi_loc[row]; 186 pa = pa_loc + pi_loc[row]; 187 188 /* perform dense axpy */ 189 valtmp = ada[j]; 190 for (k=0; k<pnz; k++){ 191 apa[pj[k]] += valtmp*pa[k]; 192 } 193 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 194 } 195 196 /* off-diagonal portion of A */ 197 anz = aoi[i+1] - aoi[i]; 198 aoj = ao->j + aoi[i]; 199 aoa = ao->a + aoi[i]; 200 for (j=0; j<anz; j++) { 201 row = aoj[j]; 202 pnz = pi_oth[row+1] - pi_oth[row]; 203 pj = pj_oth + pi_oth[row]; 204 pa = pa_oth + pi_oth[row]; 205 206 /* perform dense axpy */ 207 valtmp = aoa[j]; 208 for (k=0; k<pnz; k++){ 209 apa[pj[k]] += valtmp*pa[k]; 210 } 211 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 212 } 213 214 /* set values in C */ 215 row = rstart + i; 216 api = ptap->api; 217 apj = ptap->apj; 218 apJ = apj + api[i]; 219 cdnz = cd->i[i+1] - cd->i[i]; 220 conz = co->i[i+1] - co->i[i]; 221 222 /* 1st off-diagoanl part of C */ 223 ca = coa + co->i[i]; 224 k = 0; 225 for (k0=0; k0<conz; k0++){ 226 if (apJ[k] >= cstart) break; 227 ca[k0] = apa[apJ[k]]; 228 apa[apJ[k]] = 0.0; 229 k++; 230 } 231 232 /* diagonal part of C */ 233 ca = cda + cd->i[i]; 234 for (k1=0; k1<cdnz; k1++){ 235 ca[k1] = apa[apJ[k]]; 236 apa[apJ[k]] = 0.0; 237 k++; 238 } 239 240 /* 2nd off-diagoanl part of C */ 241 ca = coa + co->i[i]; 242 for (; k0<conz; k0++){ 243 ca[k0] = apa[apJ[k]]; 244 apa[apJ[k]] = 0.0; 245 k++; 246 } 247 } 248 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 249 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 250 #if defined(DEBUG_MATMATMULT) 251 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank); 252 #endif 253 PetscFunctionReturn(0); 254 } 255 256 #undef __FUNCT__ 257 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 258 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 259 { 260 PetscErrorCode ierr; 261 MPI_Comm comm=((PetscObject)A)->comm; 262 Mat Cmpi; 263 Mat_PtAPMPI *ptap; 264 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 265 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 266 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 267 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 268 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 269 PetscInt nlnk,*lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 270 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 271 PetscBT lnkbt; 272 PetscScalar *apa; 273 PetscReal afill; 274 PetscBool matmatmult_old=PETSC_FALSE; 275 #if defined(DEBUG_MATMATMULT) 276 PetscMPIInt rank=a->rank; 277 #endif 278 279 PetscFunctionBegin; 280 if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){ 281 SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend); 282 } 283 284 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 285 ierr = PetscOptionsBool("-matmatmult_32","Use PETSc-3.2","",matmatmult_old,&matmatmult_old,PETSC_NULL);CHKERRQ(ierr); 286 if (matmatmult_old){ 287 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_32(A,P,fill,C);CHKERRQ(ierr); 288 PetscFunctionReturn(0); 289 } 290 ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",matmatmult_old,&matmatmult_old,PETSC_NULL);CHKERRQ(ierr); 291 if (matmatmult_old){ 292 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr); 293 PetscFunctionReturn(0); 294 } 295 ierr = PetscOptionsEnd();CHKERRQ(ierr); 296 297 #if defined(DEBUG_MATMATMULT) 298 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank); 299 #endif 300 301 /* create struct Mat_PtAPMPI and attached it to C later */ 302 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 303 304 /* malloc apa to store dense row A[i,:]*P */ 305 ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 306 ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr); 307 ptap->apa = apa; 308 #if defined(DEBUG_MATMATMULT) 309 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN); 310 //if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Annz %D\n",rank,ad->i[am]+ao->i[am]); 311 #endif 312 313 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 314 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 315 #if defined(DEBUG_MATMATMULT) 316 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 317 #endif 318 /* get P_loc by taking all local rows of P */ 319 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 320 321 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 322 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 323 pi_loc = p_loc->i; pj_loc = p_loc->j; 324 pi_oth = p_oth->i; pj_oth = p_oth->j; 325 326 #if defined(DEBUG_MATMATMULT) 327 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done, Annz %D * P_locnnz %D = %D...\n",rank,ad->i[am]+ao->i[am],p_loc->rmax,(ad->i[am]+ao->i[am])*p_loc->rmax); 328 #endif 329 330 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 331 /*-------------------------------------------------------------------*/ 332 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 333 ptap->api = api; 334 api[0] = 0; 335 336 /* create and initialize a linked list */ 337 nlnk = pN+1; 338 ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 339 340 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 341 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 342 current_space = free_space; 343 344 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 345 for (i=0; i<am; i++) { 346 apnz = 0; 347 /* diagonal portion of A */ 348 nzi = adi[i+1] - adi[i]; 349 for (j=0; j<nzi; j++){ 350 row = *adj++; 351 pnz = pi_loc[row+1] - pi_loc[row]; 352 Jptr = pj_loc + pi_loc[row]; 353 /* add non-zero cols of P into the sorted linked list lnk */ 354 ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 355 apnz += nlnk; 356 } 357 /* off-diagonal portion of A */ 358 nzi = aoi[i+1] - aoi[i]; 359 for (j=0; j<nzi; j++){ 360 row = *aoj++; 361 pnz = pi_oth[row+1] - pi_oth[row]; 362 Jptr = pj_oth + pi_oth[row]; 363 ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 364 apnz += nlnk; 365 } 366 367 api[i+1] = api[i] + apnz; 368 369 /* if free space is not available, double the total space in the list */ 370 if (current_space->local_remaining<apnz) { 371 ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 372 nspacedouble++; 373 } 374 375 /* Copy data into free space, then initialize lnk */ 376 ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 377 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 378 current_space->array += apnz; 379 current_space->local_used += apnz; 380 current_space->local_remaining -= apnz; 381 } 382 383 /* Allocate space for apj, initialize apj, and */ 384 /* destroy list of free space and other temporary array(s) */ 385 ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 386 apj = ptap->apj; 387 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 388 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 389 #if defined(DEBUG_MATMATMULT) 390 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done...\n",rank); 391 #endif 392 393 /* create and assemble symbolic parallel matrix Cmpi */ 394 /*----------------------------------------------------*/ 395 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 396 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 397 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 398 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 399 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 400 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 401 for (i=0; i<am; i++){ 402 row = i + rstart; 403 apnz = api[i+1] - api[i]; 404 ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 405 apj += apnz; 406 } 407 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 408 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 409 410 ptap->destroy = Cmpi->ops->destroy; 411 ptap->duplicate = Cmpi->ops->duplicate; 412 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 413 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 414 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 415 416 /* attach the supporting struct to Cmpi for reuse */ 417 c = (Mat_MPIAIJ*)Cmpi->data; 418 c->ptap = ptap; 419 420 *C = Cmpi; 421 422 /* set MatInfo */ 423 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 424 if (afill < 1.0) afill = 1.0; 425 Cmpi->info.mallocs = nspacedouble; 426 Cmpi->info.fill_ratio_given = fill; 427 Cmpi->info.fill_ratio_needed = afill; 428 429 #if defined(PETSC_USE_INFO) 430 if (api[am]) { 431 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 432 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 433 } else { 434 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 435 } 436 #endif 437 #if defined(DEBUG_MATMATMULT) 438 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank); 439 #endif 440 PetscFunctionReturn(0); 441 } 442 443 /* implementation used in PETSc-3.2 */ 444 /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 445 #undef __FUNCT__ 446 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 447 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat B,Mat C) 448 { 449 PetscErrorCode ierr; 450 Mat *seq; 451 Mat_MatMatMultMPI *mult; 452 PetscContainer container; 453 454 PetscFunctionBegin; 455 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 456 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 457 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 458 459 if (mult->skipNumeric){ /* first numeric product is done during symbolic product */ 460 mult->skipNumeric = PETSC_FALSE; 461 PetscFunctionReturn(0); 462 } 463 #if defined(DEBUG_MATMATMULT) 464 PetscMPIInt rank; 465 MPI_Comm comm = ((PetscObject)C)->comm; 466 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 467 ierr = MPI_Barrier(comm);CHKERRQ(ierr); 468 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 469 #endif 470 471 seq = &mult->B_seq; 472 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 473 mult->B_seq = *seq; 474 475 seq = &mult->A_loc; 476 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 477 mult->A_loc = *seq; 478 479 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 480 481 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 482 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 483 #if defined(DEBUG_MATMATMULT) 484 ierr = MPI_Barrier(comm);CHKERRQ(ierr); 485 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 486 #endif 487 PetscFunctionReturn(0); 488 } 489 490 /* numeric product is computed as well */ 491 #undef __FUNCT__ 492 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 493 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 494 { 495 PetscErrorCode ierr; 496 Mat_MatMatMultMPI *mult; 497 PetscContainer container; 498 Mat AB,*seq; 499 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 500 PetscInt *idx,i,start,ncols,nzA,nzB,*cmap,imark; 501 #if defined(DEBUG_MATMATMULT) 502 MPI_Comm comm = ((PetscObject)A)->comm; 503 PetscMPIInt rank=a->rank; 504 #endif 505 506 PetscFunctionBegin; 507 #if defined(DEBUG_MATMATMULT) 508 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 509 #endif 510 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 511 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); 512 } 513 514 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 515 516 /* get isrowb: nonzero col of A */ 517 start = A->cmap->rstart; 518 cmap = a->garray; 519 nzA = a->A->cmap->n; 520 nzB = a->B->cmap->n; 521 ierr = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr); 522 ncols = 0; 523 for (i=0; i<nzB; i++) { /* row < local row index */ 524 if (cmap[i] < start) idx[ncols++] = cmap[i]; 525 else break; 526 } 527 imark = i; 528 for (i=0; i<nzA; i++) idx[ncols++] = start + i; /* local rows */ 529 for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */ 530 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr); 531 ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr); 532 533 /* get isrowa: all local rows of A */ 534 ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr); 535 536 /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */ 537 /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 538 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 539 mult->B_seq = *seq; 540 ierr = PetscFree(seq);CHKERRQ(ierr); 541 #if defined(DEBUG_MATMATMULT) 542 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] B_seq is done...\n",rank); 543 #endif 544 545 /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 546 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 547 mult->A_loc = *seq; 548 ierr = PetscFree(seq);CHKERRQ(ierr); 549 #if defined(DEBUG_MATMATMULT) 550 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A_loc is done...\n",rank); 551 #endif 552 553 /* compute C_seq = A_seq * B_seq */ 554 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr); 555 #if defined(DEBUG_MATMATMULT) 556 ierr = MPI_Barrier(comm);CHKERRQ(ierr); 557 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Symbolic is done...\n",rank); 558 #endif 559 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 560 #if defined(DEBUG_MATMATMULT) 561 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Numeric is done...\n",rank); 562 #endif 563 564 /* create mpi matrix C by concatinating C_seq */ 565 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 566 ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr); 567 ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr); 568 #if defined(DEBUG_MATMATMULT) 569 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank); 570 #endif 571 572 /* attach the supporting struct to C for reuse of symbolic C */ 573 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 574 ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 575 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 576 ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 577 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 578 mult->skipNumeric = PETSC_TRUE; /* a numeric product is done here */ 579 mult->destroy = AB->ops->destroy; 580 mult->duplicate = AB->ops->duplicate; 581 AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 582 AB->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_32; 583 AB->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_32; 584 AB->ops->matmult = MatMatMult_MPIAIJ_MPIAIJ; 585 *C = AB; 586 #if defined(DEBUG_MATMATMULT) 587 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 588 #endif 589 PetscFunctionReturn(0); 590 } 591 592 #undef __FUNCT__ 593 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 594 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 595 { 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 if (scall == MAT_INITIAL_MATRIX){ 600 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 601 } 602 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 603 PetscFunctionReturn(0); 604 } 605 606 typedef struct { 607 Mat workB; 608 PetscScalar *rvalues,*svalues; 609 MPI_Request *rwaits,*swaits; 610 } MPIAIJ_MPIDense; 611 612 #undef __FUNCT__ 613 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 614 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 615 { 616 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 617 PetscErrorCode ierr; 618 619 PetscFunctionBegin; 620 ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 621 ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 622 ierr = PetscFree(contents);CHKERRQ(ierr); 623 PetscFunctionReturn(0); 624 } 625 626 #undef __FUNCT__ 627 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 628 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 629 { 630 PetscErrorCode ierr; 631 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 632 PetscInt nz = aij->B->cmap->n; 633 PetscContainer container; 634 MPIAIJ_MPIDense *contents; 635 VecScatter ctx = aij->Mvctx; 636 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 637 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 638 PetscInt m=A->rmap->n,n=B->cmap->n; 639 640 PetscFunctionBegin; 641 ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 642 ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 643 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 644 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 645 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 646 (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 647 648 ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 649 /* Create work matrix used to store off processor rows of B needed for local product */ 650 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 651 /* Create work arrays needed */ 652 ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 653 B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 654 from->n,MPI_Request,&contents->rwaits, 655 to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 656 657 ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 658 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 659 ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 660 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 661 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 662 PetscFunctionReturn(0); 663 } 664 665 #undef __FUNCT__ 666 #define __FUNCT__ "MatMPIDenseScatter" 667 /* 668 Performs an efficient scatter on the rows of B needed by this process; this is 669 a modification of the VecScatterBegin_() routines. 670 */ 671 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 672 { 673 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 674 PetscErrorCode ierr; 675 PetscScalar *b,*w,*svalues,*rvalues; 676 VecScatter ctx = aij->Mvctx; 677 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 678 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 679 PetscInt i,j,k; 680 PetscInt *sindices,*sstarts,*rindices,*rstarts; 681 PetscMPIInt *sprocs,*rprocs,nrecvs; 682 MPI_Request *swaits,*rwaits; 683 MPI_Comm comm = ((PetscObject)A)->comm; 684 PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 685 MPI_Status status; 686 MPIAIJ_MPIDense *contents; 687 PetscContainer container; 688 Mat workB; 689 690 PetscFunctionBegin; 691 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 692 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 693 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 694 695 workB = *outworkB = contents->workB; 696 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); 697 sindices = to->indices; 698 sstarts = to->starts; 699 sprocs = to->procs; 700 swaits = contents->swaits; 701 svalues = contents->svalues; 702 703 rindices = from->indices; 704 rstarts = from->starts; 705 rprocs = from->procs; 706 rwaits = contents->rwaits; 707 rvalues = contents->rvalues; 708 709 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 710 ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 711 712 for (i=0; i<from->n; i++) { 713 ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 714 } 715 716 for (i=0; i<to->n; i++) { 717 /* pack a message at a time */ 718 CHKMEMQ; 719 for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 720 for (k=0; k<ncols; k++) { 721 svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 722 } 723 } 724 CHKMEMQ; 725 ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 726 } 727 728 nrecvs = from->n; 729 while (nrecvs) { 730 ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 731 nrecvs--; 732 /* unpack a message at a time */ 733 CHKMEMQ; 734 for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 735 for (k=0; k<ncols; k++) { 736 w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 737 } 738 } 739 CHKMEMQ; 740 } 741 if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 742 743 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 744 ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 745 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 746 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 747 PetscFunctionReturn(0); 748 } 749 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 750 751 #undef __FUNCT__ 752 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 753 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 754 { 755 PetscErrorCode ierr; 756 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 757 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 758 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 759 Mat workB; 760 761 PetscFunctionBegin; 762 763 /* diagonal block of A times all local rows of B*/ 764 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 765 766 /* get off processor parts of B needed to complete the product */ 767 ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 768 769 /* off-diagonal block of A times nonlocal rows of B */ 770 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 771 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 772 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 773 PetscFunctionReturn(0); 774 } 775 /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 776 #undef __FUNCT__ 777 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_32" 778 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_32(Mat A,Mat B,Mat C) 779 { 780 PetscErrorCode ierr; 781 Mat *seq; 782 Mat_MatMatMultMPI *mult; 783 PetscContainer container; 784 785 PetscFunctionBegin; 786 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 787 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 788 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 789 if (mult->skipNumeric){ /* first numeric product is done during symbolic product */ 790 mult->skipNumeric = PETSC_FALSE; 791 PetscFunctionReturn(0); 792 } 793 #if defined(DEBUG_MATMATMULT) 794 PetscMPIInt rank; 795 ierr = MPI_Comm_rank(((PetscObject)C)->comm,&rank);CHKERRQ(ierr); 796 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_32()...\n",rank); 797 #endif 798 799 seq = &mult->B_seq; 800 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 801 mult->B_seq = *seq; 802 803 seq = &mult->A_loc; 804 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 805 mult->A_loc = *seq; 806 807 ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 808 809 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 810 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 811 PetscFunctionReturn(0); 812 } 813 814 #undef __FUNCT__ 815 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_32" 816 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_32(Mat A,Mat B,PetscReal fill,Mat *C) 817 { 818 PetscErrorCode ierr; 819 PetscInt start,end; 820 Mat_MatMatMultMPI *mult; 821 PetscContainer container; 822 823 PetscFunctionBegin; 824 #if defined(DEBUG_MATMATMULT) 825 PetscMPIInt rank; 826 ierr = MPI_Comm_rank(((PetscObject)A)->comm,&rank);CHKERRQ(ierr); 827 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_32()...\n",rank); 828 #endif 829 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 830 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); 831 } 832 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 833 834 /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 835 ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->B_seq);CHKERRQ(ierr); 836 837 /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 838 start = A->rmap->rstart; end = A->rmap->rend; 839 ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr); 840 ierr = MatMPIAIJGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr); 841 842 /* compute C_seq = A_seq * B_seq */ 843 ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr); 844 845 /* create mpi matrix C by concatinating C_seq */ 846 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 847 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 848 849 /* attach the supporting struct to C for reuse of symbolic C */ 850 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 851 ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 852 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 853 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 854 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 855 mult->skipNumeric = PETSC_TRUE; /* a numeric product is done here */ 856 mult->destroy = (*C)->ops->destroy; 857 mult->duplicate = (*C)->ops->duplicate; 858 (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_32; 859 (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_32; 860 (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_32; 861 PetscFunctionReturn(0); 862 } 863