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 ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_32",&matmatmult_old,PETSC_NULL);CHKERRQ(ierr); 284 if (matmatmult_old){ 285 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_32(A,P,fill,C);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_scalable",&matmatmult_old,PETSC_NULL);CHKERRQ(ierr); 289 if (matmatmult_old){ 290 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 #if defined(DEBUG_MATMATMULT) 295 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank); 296 #endif 297 298 /* create struct Mat_PtAPMPI and attached it to C later */ 299 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 300 301 /* malloc apa to store dense row A[i,:]*P */ 302 ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 303 ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr); 304 ptap->apa = apa; 305 #if defined(DEBUG_MATMATMULT) 306 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN); 307 //if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Annz %D\n",rank,ad->i[am]+ao->i[am]); 308 #endif 309 310 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 311 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 312 #if defined(DEBUG_MATMATMULT) 313 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 314 #endif 315 /* get P_loc by taking all local rows of P */ 316 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 317 318 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 319 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 320 pi_loc = p_loc->i; pj_loc = p_loc->j; 321 pi_oth = p_oth->i; pj_oth = p_oth->j; 322 323 #if defined(DEBUG_MATMATMULT) 324 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); 325 #endif 326 327 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 328 /*-------------------------------------------------------------------*/ 329 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 330 ptap->api = api; 331 api[0] = 0; 332 333 /* create and initialize a linked list */ 334 nlnk = pN+1; 335 ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 336 337 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 338 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 339 current_space = free_space; 340 341 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 342 for (i=0; i<am; i++) { 343 apnz = 0; 344 /* diagonal portion of A */ 345 nzi = adi[i+1] - adi[i]; 346 for (j=0; j<nzi; j++){ 347 row = *adj++; 348 pnz = pi_loc[row+1] - pi_loc[row]; 349 Jptr = pj_loc + pi_loc[row]; 350 /* add non-zero cols of P into the sorted linked list lnk */ 351 ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 352 apnz += nlnk; 353 } 354 /* off-diagonal portion of A */ 355 nzi = aoi[i+1] - aoi[i]; 356 for (j=0; j<nzi; j++){ 357 row = *aoj++; 358 pnz = pi_oth[row+1] - pi_oth[row]; 359 Jptr = pj_oth + pi_oth[row]; 360 ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 361 apnz += nlnk; 362 } 363 364 api[i+1] = api[i] + apnz; 365 366 /* if free space is not available, double the total space in the list */ 367 if (current_space->local_remaining<apnz) { 368 ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 369 nspacedouble++; 370 } 371 372 /* Copy data into free space, then initialize lnk */ 373 ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 374 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 375 current_space->array += apnz; 376 current_space->local_used += apnz; 377 current_space->local_remaining -= apnz; 378 } 379 380 /* Allocate space for apj, initialize apj, and */ 381 /* destroy list of free space and other temporary array(s) */ 382 ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 383 apj = ptap->apj; 384 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 385 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 386 #if defined(DEBUG_MATMATMULT) 387 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done...\n",rank); 388 #endif 389 390 /* create and assemble symbolic parallel matrix Cmpi */ 391 /*----------------------------------------------------*/ 392 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 393 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 394 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 395 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 396 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 397 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 398 for (i=0; i<am; i++){ 399 row = i + rstart; 400 apnz = api[i+1] - api[i]; 401 ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 402 apj += apnz; 403 } 404 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 405 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 406 407 ptap->destroy = Cmpi->ops->destroy; 408 ptap->duplicate = Cmpi->ops->duplicate; 409 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 410 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 411 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 412 413 /* attach the supporting struct to Cmpi for reuse */ 414 c = (Mat_MPIAIJ*)Cmpi->data; 415 c->ptap = ptap; 416 417 *C = Cmpi; 418 419 /* set MatInfo */ 420 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 421 if (afill < 1.0) afill = 1.0; 422 Cmpi->info.mallocs = nspacedouble; 423 Cmpi->info.fill_ratio_given = fill; 424 Cmpi->info.fill_ratio_needed = afill; 425 426 #if defined(PETSC_USE_INFO) 427 if (api[am]) { 428 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 429 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 430 } else { 431 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 432 } 433 #endif 434 #if defined(DEBUG_MATMATMULT) 435 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank); 436 #endif 437 PetscFunctionReturn(0); 438 } 439 440 /* implementation used in PETSc-3.2 */ 441 /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 442 #undef __FUNCT__ 443 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 444 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat B,Mat C) 445 { 446 PetscErrorCode ierr; 447 Mat *seq; 448 Mat_MatMatMultMPI *mult; 449 PetscContainer container; 450 451 PetscFunctionBegin; 452 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 453 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 454 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 455 456 if (mult->skipNumeric){ /* first numeric product is done during symbolic product */ 457 mult->skipNumeric = PETSC_FALSE; 458 PetscFunctionReturn(0); 459 } 460 #if defined(DEBUG_MATMATMULT) 461 PetscMPIInt rank; 462 MPI_Comm comm = ((PetscObject)C)->comm; 463 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 464 ierr = MPI_Barrier(comm);CHKERRQ(ierr); 465 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 466 #endif 467 468 seq = &mult->B_seq; 469 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 470 mult->B_seq = *seq; 471 472 seq = &mult->A_loc; 473 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 474 mult->A_loc = *seq; 475 476 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 477 478 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 479 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 480 #if defined(DEBUG_MATMATMULT) 481 ierr = MPI_Barrier(comm);CHKERRQ(ierr); 482 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 483 #endif 484 PetscFunctionReturn(0); 485 } 486 487 /* numeric product is computed as well */ 488 #undef __FUNCT__ 489 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 490 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 491 { 492 PetscErrorCode ierr; 493 Mat_MatMatMultMPI *mult; 494 PetscContainer container; 495 Mat AB,*seq; 496 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 497 PetscInt *idx,i,start,ncols,nzA,nzB,*cmap,imark; 498 #if defined(DEBUG_MATMATMULT) 499 MPI_Comm comm = ((PetscObject)A)->comm; 500 PetscMPIInt rank=a->rank; 501 #endif 502 503 PetscFunctionBegin; 504 #if defined(DEBUG_MATMATMULT) 505 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 506 #endif 507 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 508 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); 509 } 510 511 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 512 513 /* get isrowb: nonzero col of A */ 514 start = A->cmap->rstart; 515 cmap = a->garray; 516 nzA = a->A->cmap->n; 517 nzB = a->B->cmap->n; 518 ierr = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr); 519 ncols = 0; 520 for (i=0; i<nzB; i++) { /* row < local row index */ 521 if (cmap[i] < start) idx[ncols++] = cmap[i]; 522 else break; 523 } 524 imark = i; 525 for (i=0; i<nzA; i++) idx[ncols++] = start + i; /* local rows */ 526 for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */ 527 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr); 528 ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr); 529 530 /* get isrowa: all local rows of A */ 531 ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr); 532 533 /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */ 534 /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 535 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 536 mult->B_seq = *seq; 537 ierr = PetscFree(seq);CHKERRQ(ierr); 538 #if defined(DEBUG_MATMATMULT) 539 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] B_seq is done...\n",rank); 540 #endif 541 542 /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 543 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 544 mult->A_loc = *seq; 545 ierr = PetscFree(seq);CHKERRQ(ierr); 546 #if defined(DEBUG_MATMATMULT) 547 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A_loc is done...\n",rank); 548 #endif 549 550 /* compute C_seq = A_seq * B_seq */ 551 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr); 552 #if defined(DEBUG_MATMATMULT) 553 ierr = MPI_Barrier(comm);CHKERRQ(ierr); 554 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Symbolic is done...\n",rank); 555 #endif 556 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 557 #if defined(DEBUG_MATMATMULT) 558 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Numeric is done...\n",rank); 559 #endif 560 561 /* create mpi matrix C by concatinating C_seq */ 562 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 563 ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr); 564 ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr); 565 #if defined(DEBUG_MATMATMULT) 566 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank); 567 #endif 568 569 /* attach the supporting struct to C for reuse of symbolic C */ 570 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 571 ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 572 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 573 ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 574 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 575 mult->skipNumeric = PETSC_TRUE; /* a numeric product is done here */ 576 mult->destroy = AB->ops->destroy; 577 mult->duplicate = AB->ops->duplicate; 578 AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 579 AB->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_32; 580 AB->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_32; 581 AB->ops->matmult = MatMatMult_MPIAIJ_MPIAIJ; 582 *C = AB; 583 #if defined(DEBUG_MATMATMULT) 584 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 585 #endif 586 PetscFunctionReturn(0); 587 } 588 589 #undef __FUNCT__ 590 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 591 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 592 { 593 PetscErrorCode ierr; 594 595 PetscFunctionBegin; 596 if (scall == MAT_INITIAL_MATRIX){ 597 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 598 } 599 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 600 PetscFunctionReturn(0); 601 } 602 603 typedef struct { 604 Mat workB; 605 PetscScalar *rvalues,*svalues; 606 MPI_Request *rwaits,*swaits; 607 } MPIAIJ_MPIDense; 608 609 #undef __FUNCT__ 610 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 611 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 612 { 613 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 614 PetscErrorCode ierr; 615 616 PetscFunctionBegin; 617 ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 618 ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 619 ierr = PetscFree(contents);CHKERRQ(ierr); 620 PetscFunctionReturn(0); 621 } 622 623 #undef __FUNCT__ 624 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 625 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 626 { 627 PetscErrorCode ierr; 628 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 629 PetscInt nz = aij->B->cmap->n; 630 PetscContainer container; 631 MPIAIJ_MPIDense *contents; 632 VecScatter ctx = aij->Mvctx; 633 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 634 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 635 PetscInt m=A->rmap->n,n=B->cmap->n; 636 637 PetscFunctionBegin; 638 ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 639 ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 640 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 641 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 642 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 643 (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 644 645 ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 646 /* Create work matrix used to store off processor rows of B needed for local product */ 647 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 648 /* Create work arrays needed */ 649 ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 650 B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 651 from->n,MPI_Request,&contents->rwaits, 652 to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 653 654 ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 655 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 656 ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 657 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 658 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 659 PetscFunctionReturn(0); 660 } 661 662 #undef __FUNCT__ 663 #define __FUNCT__ "MatMPIDenseScatter" 664 /* 665 Performs an efficient scatter on the rows of B needed by this process; this is 666 a modification of the VecScatterBegin_() routines. 667 */ 668 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 669 { 670 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 671 PetscErrorCode ierr; 672 PetscScalar *b,*w,*svalues,*rvalues; 673 VecScatter ctx = aij->Mvctx; 674 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 675 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 676 PetscInt i,j,k; 677 PetscInt *sindices,*sstarts,*rindices,*rstarts; 678 PetscMPIInt *sprocs,*rprocs,nrecvs; 679 MPI_Request *swaits,*rwaits; 680 MPI_Comm comm = ((PetscObject)A)->comm; 681 PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 682 MPI_Status status; 683 MPIAIJ_MPIDense *contents; 684 PetscContainer container; 685 Mat workB; 686 687 PetscFunctionBegin; 688 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 689 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 690 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 691 692 workB = *outworkB = contents->workB; 693 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); 694 sindices = to->indices; 695 sstarts = to->starts; 696 sprocs = to->procs; 697 swaits = contents->swaits; 698 svalues = contents->svalues; 699 700 rindices = from->indices; 701 rstarts = from->starts; 702 rprocs = from->procs; 703 rwaits = contents->rwaits; 704 rvalues = contents->rvalues; 705 706 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 707 ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 708 709 for (i=0; i<from->n; i++) { 710 ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 711 } 712 713 for (i=0; i<to->n; i++) { 714 /* pack a message at a time */ 715 CHKMEMQ; 716 for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 717 for (k=0; k<ncols; k++) { 718 svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 719 } 720 } 721 CHKMEMQ; 722 ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 723 } 724 725 nrecvs = from->n; 726 while (nrecvs) { 727 ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 728 nrecvs--; 729 /* unpack a message at a time */ 730 CHKMEMQ; 731 for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 732 for (k=0; k<ncols; k++) { 733 w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 734 } 735 } 736 CHKMEMQ; 737 } 738 if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 739 740 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 741 ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 742 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 743 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 744 PetscFunctionReturn(0); 745 } 746 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 747 748 #undef __FUNCT__ 749 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 750 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 751 { 752 PetscErrorCode ierr; 753 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 754 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 755 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 756 Mat workB; 757 758 PetscFunctionBegin; 759 760 /* diagonal block of A times all local rows of B*/ 761 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 762 763 /* get off processor parts of B needed to complete the product */ 764 ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 765 766 /* off-diagonal block of A times nonlocal rows of B */ 767 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 768 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 769 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 770 PetscFunctionReturn(0); 771 } 772 /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 773 #undef __FUNCT__ 774 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_32" 775 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_32(Mat A,Mat B,Mat C) 776 { 777 PetscErrorCode ierr; 778 Mat *seq; 779 Mat_MatMatMultMPI *mult; 780 PetscContainer container; 781 782 PetscFunctionBegin; 783 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 784 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 785 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 786 if (mult->skipNumeric){ /* first numeric product is done during symbolic product */ 787 mult->skipNumeric = PETSC_FALSE; 788 PetscFunctionReturn(0); 789 } 790 #if defined(DEBUG_MATMATMULT) 791 PetscMPIInt rank; 792 ierr = MPI_Comm_rank(((PetscObject)C)->comm,&rank);CHKERRQ(ierr); 793 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_32()...\n",rank); 794 #endif 795 796 seq = &mult->B_seq; 797 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 798 mult->B_seq = *seq; 799 800 seq = &mult->A_loc; 801 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 802 mult->A_loc = *seq; 803 804 ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 805 806 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 807 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 808 PetscFunctionReturn(0); 809 } 810 811 #undef __FUNCT__ 812 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_32" 813 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_32(Mat A,Mat B,PetscReal fill,Mat *C) 814 { 815 PetscErrorCode ierr; 816 PetscInt start,end; 817 Mat_MatMatMultMPI *mult; 818 PetscContainer container; 819 820 PetscFunctionBegin; 821 #if defined(DEBUG_MATMATMULT) 822 PetscMPIInt rank; 823 ierr = MPI_Comm_rank(((PetscObject)A)->comm,&rank);CHKERRQ(ierr); 824 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_32()...\n",rank); 825 #endif 826 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 827 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); 828 } 829 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 830 831 /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 832 ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->B_seq);CHKERRQ(ierr); 833 834 /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 835 start = A->rmap->rstart; end = A->rmap->rend; 836 ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr); 837 ierr = MatMPIAIJGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr); 838 839 /* compute C_seq = A_seq * B_seq */ 840 ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr); 841 842 /* create mpi matrix C by concatinating C_seq */ 843 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 844 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 845 846 /* attach the supporting struct to C for reuse of symbolic C */ 847 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 848 ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 849 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 850 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 851 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 852 mult->skipNumeric = PETSC_TRUE; /* a numeric product is done here */ 853 mult->destroy = (*C)->ops->destroy; 854 mult->duplicate = (*C)->ops->duplicate; 855 (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_32; 856 (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_32; 857 (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_32; 858 PetscFunctionReturn(0); 859 } 860