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 #define DEBUG_MATTrMatMult 16 */ 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 20 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 21 { 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 if (scall == MAT_INITIAL_MATRIX){ 26 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 27 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 28 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 29 } 30 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 31 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 32 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI" 38 PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr) 39 { 40 PetscErrorCode ierr; 41 Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)ptr; 42 43 PetscFunctionBegin; 44 ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr); 45 ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr); 46 ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr); 47 ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr); 48 ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr); 49 ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr); 50 ierr = PetscFree(mult);CHKERRQ(ierr); 51 PetscFunctionReturn(0); 52 } 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 56 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 57 { 58 PetscErrorCode ierr; 59 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 60 Mat_PtAPMPI *ptap=a->ptap; 61 62 PetscFunctionBegin; 63 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 64 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 65 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 66 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 67 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 68 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 69 ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 70 ierr = ptap->destroy(A);CHKERRQ(ierr); 71 ierr = PetscFree(ptap);CHKERRQ(ierr); 72 PetscFunctionReturn(0); 73 } 74 75 #undef __FUNCT__ 76 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32" 77 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A) 78 { 79 PetscErrorCode ierr; 80 PetscContainer container; 81 Mat_MatMatMultMPI *mult=PETSC_NULL; 82 83 PetscFunctionBegin; 84 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 85 if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 86 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 87 A->ops->destroy = mult->destroy; 88 A->ops->duplicate = mult->duplicate; 89 if (A->ops->destroy) { 90 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 91 } 92 ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32" 98 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M) 99 { 100 PetscErrorCode ierr; 101 Mat_MatMatMultMPI *mult; 102 PetscContainer container; 103 104 PetscFunctionBegin; 105 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 106 if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 107 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 108 /* Note: the container is not duplicated, because it requires deep copying of 109 several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()). 110 These data sets are only used for repeated calling of MatMatMultNumeric(). 111 *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */ 112 ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr); 113 (*M)->ops->destroy = mult->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 114 (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */ 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 120 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 121 { 122 PetscErrorCode ierr; 123 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 124 Mat_PtAPMPI *ptap=a->ptap; 125 126 PetscFunctionBegin; 127 ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 128 (*M)->ops->destroy = ptap->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */ 129 (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */ 130 PetscFunctionReturn(0); 131 } 132 133 #undef __FUNCT__ 134 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 135 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 136 { 137 PetscErrorCode ierr; 138 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 139 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 140 Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 141 PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 142 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 143 Mat_SeqAIJ *p_loc,*p_oth; 144 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 145 PetscScalar *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca; 146 PetscInt cm=C->rmap->n,anz,pnz; 147 Mat_PtAPMPI *ptap=c->ptap; 148 PetscInt *api,*apj,*apJ,i,j,k,row; 149 PetscInt cstart=C->cmap->rstart; 150 PetscInt cdnz,conz,k0,k1; 151 #if defined(DEBUG_MATMATMULT) 152 PetscMPIInt rank=a->rank; 153 #endif 154 155 PetscFunctionBegin; 156 #if defined(DEBUG_MATMATMULT) 157 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank); 158 #endif 159 160 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 161 /*-----------------------------------------------------*/ 162 /* update numerical values of P_oth and P_loc */ 163 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 164 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 165 #if defined(DEBUG_MATMATMULT) 166 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 167 #endif 168 169 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 170 /*----------------------------------------------------------*/ 171 /* get data from symbolic products */ 172 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 173 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 174 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 175 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 176 177 /* get apa for storing dense row A[i,:]*P */ 178 apa = ptap->apa; 179 180 api = ptap->api; 181 apj = ptap->apj; 182 for (i=0; i<cm; i++) { 183 /* diagonal portion of A */ 184 anz = adi[i+1] - adi[i]; 185 adj = ad->j + adi[i]; 186 ada = ad->a + adi[i]; 187 for (j=0; j<anz; j++) { 188 row = adj[j]; 189 pnz = pi_loc[row+1] - pi_loc[row]; 190 pj = pj_loc + pi_loc[row]; 191 pa = pa_loc + pi_loc[row]; 192 193 /* perform dense axpy */ 194 valtmp = ada[j]; 195 for (k=0; k<pnz; k++){ 196 apa[pj[k]] += valtmp*pa[k]; 197 } 198 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 199 } 200 201 /* off-diagonal portion of A */ 202 anz = aoi[i+1] - aoi[i]; 203 aoj = ao->j + aoi[i]; 204 aoa = ao->a + aoi[i]; 205 for (j=0; j<anz; j++) { 206 row = aoj[j]; 207 pnz = pi_oth[row+1] - pi_oth[row]; 208 pj = pj_oth + pi_oth[row]; 209 pa = pa_oth + pi_oth[row]; 210 211 /* perform dense axpy */ 212 valtmp = aoa[j]; 213 for (k=0; k<pnz; k++){ 214 apa[pj[k]] += valtmp*pa[k]; 215 } 216 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 217 } 218 219 /* set values in C */ 220 apJ = apj + api[i]; 221 cdnz = cd->i[i+1] - cd->i[i]; 222 conz = co->i[i+1] - co->i[i]; 223 224 /* 1st off-diagoanl part of C */ 225 ca = coa + co->i[i]; 226 k = 0; 227 for (k0=0; k0<conz; k0++){ 228 if (apJ[k] >= cstart) break; 229 ca[k0] = apa[apJ[k]]; 230 apa[apJ[k]] = 0.0; 231 k++; 232 } 233 234 /* diagonal part of C */ 235 ca = cda + cd->i[i]; 236 for (k1=0; k1<cdnz; k1++){ 237 ca[k1] = apa[apJ[k]]; 238 apa[apJ[k]] = 0.0; 239 k++; 240 } 241 242 /* 2nd off-diagoanl part of C */ 243 ca = coa + co->i[i]; 244 for (; k0<conz; k0++){ 245 ca[k0] = apa[apJ[k]]; 246 apa[apJ[k]] = 0.0; 247 k++; 248 } 249 } 250 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 251 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 252 PetscFunctionReturn(0); 253 } 254 255 #undef __FUNCT__ 256 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 257 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 258 { 259 PetscErrorCode ierr; 260 MPI_Comm comm=((PetscObject)A)->comm; 261 Mat Cmpi; 262 Mat_PtAPMPI *ptap; 263 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 264 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 265 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 266 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 267 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 268 PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 269 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 270 PetscBT lnkbt; 271 PetscScalar *apa; 272 PetscReal afill; 273 PetscBool scalable=PETSC_FALSE; 274 PetscInt nlnk_max,armax,prmax; 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(comm,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 286 ierr = PetscOptionsBool("-matmatmult_32","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 287 if (scalable){ 288 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(A,P,fill,C);;CHKERRQ(ierr); 289 PetscFunctionReturn(0); 290 } 291 ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 292 if (scalable){ 293 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr); 294 PetscFunctionReturn(0); 295 } 296 ierr = PetscOptionsEnd();CHKERRQ(ierr); 297 298 #if defined(DEBUG_MATMATMULT) 299 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank); 300 #endif 301 302 /* create struct Mat_PtAPMPI and attached it to C later */ 303 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 304 305 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 306 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 307 #if defined(DEBUG_MATMATMULT) 308 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 309 #endif 310 /* get P_loc by taking all local rows of P */ 311 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 312 313 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 314 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 315 pi_loc = p_loc->i; pj_loc = p_loc->j; 316 pi_oth = p_oth->i; pj_oth = p_oth->j; 317 318 #if defined(DEBUG_MATMATMULT) 319 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); 320 #endif 321 322 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 323 /*-------------------------------------------------------------------*/ 324 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 325 ptap->api = api; 326 api[0] = 0; 327 328 /* create and initialize a linked list */ 329 armax = ad->rmax+ao->rmax; 330 prmax = PetscMax(p_loc->rmax,p_oth->rmax); 331 nlnk_max = armax*prmax; 332 if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 333 #if defined(DEBUG_MATMATMULT) 334 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] pN %d; nlnk_max %d; Armax %d+%d=%d; Prmax Max(%d,%d)=%d\n",rank,pN,nlnk_max,ad->rmax,ao->rmax,armax,p_loc->rmax,p_oth->rmax,prmax); 335 #endif 336 ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr); 337 338 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 339 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 340 current_space = free_space; 341 342 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 343 for (i=0; i<am; i++) { 344 apnz = 0; 345 /* diagonal portion of A */ 346 nzi = adi[i+1] - adi[i]; 347 for (j=0; j<nzi; j++){ 348 row = *adj++; 349 pnz = pi_loc[row+1] - pi_loc[row]; 350 Jptr = pj_loc + pi_loc[row]; 351 /* add non-zero cols of P into the sorted linked list lnk */ 352 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 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 = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 361 } 362 363 apnz = lnk[0]; 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 = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,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 /* malloc apa to store dense row A[i,:]*P */ 391 ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 392 ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr); 393 ptap->apa = apa; 394 #if defined(DEBUG_MATMATMULT) 395 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN); 396 #endif 397 398 /* create and assemble symbolic parallel matrix Cmpi */ 399 /*----------------------------------------------------*/ 400 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 401 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 402 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 403 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 404 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 405 for (i=0; i<am; i++){ 406 row = i + rstart; 407 apnz = api[i+1] - api[i]; 408 ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 409 apj += apnz; 410 } 411 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 412 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 413 414 ptap->destroy = Cmpi->ops->destroy; 415 ptap->duplicate = Cmpi->ops->duplicate; 416 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 417 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 418 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 419 420 /* attach the supporting struct to Cmpi for reuse */ 421 c = (Mat_MPIAIJ*)Cmpi->data; 422 c->ptap = ptap; 423 424 *C = Cmpi; 425 426 /* set MatInfo */ 427 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 428 if (afill < 1.0) afill = 1.0; 429 Cmpi->info.mallocs = nspacedouble; 430 Cmpi->info.fill_ratio_given = fill; 431 Cmpi->info.fill_ratio_needed = afill; 432 433 #if defined(PETSC_USE_INFO) 434 if (api[am]) { 435 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 436 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 437 } else { 438 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 439 } 440 #endif 441 PetscFunctionReturn(0); 442 } 443 444 /* implementation used in PETSc-3.2 */ 445 /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 446 #undef __FUNCT__ 447 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32" 448 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,Mat C) 449 { 450 PetscErrorCode ierr; 451 Mat *seq; 452 Mat_MatMatMultMPI *mult; 453 PetscContainer container; 454 455 PetscFunctionBegin; 456 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 457 if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 458 ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 459 460 if (mult->skipNumeric){ /* first numeric product is done during symbolic product */ 461 mult->skipNumeric = PETSC_FALSE; 462 PetscFunctionReturn(0); 463 } 464 #if defined(DEBUG_MATMATMULT) 465 PetscMPIInt rank; 466 MPI_Comm comm = ((PetscObject)C)->comm; 467 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 468 ierr = MPI_Barrier(comm);CHKERRQ(ierr); 469 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank); 470 #endif 471 472 seq = &mult->B_seq; 473 ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 474 mult->B_seq = *seq; 475 476 seq = &mult->A_loc; 477 ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 478 mult->A_loc = *seq; 479 480 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 481 482 ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 483 ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 484 PetscFunctionReturn(0); 485 } 486 487 /* numeric product is computed as well */ 488 #undef __FUNCT__ 489 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32" 490 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(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_32()...\n",rank); 506 #endif 507 if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 508 SETERRQ4(((PetscObject)A)->comm,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 concatenating 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 ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr); 566 #if defined(DEBUG_MATMATMULT) 567 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank); 568 #endif 569 570 /* attach the supporting struct to C for reuse of symbolic C */ 571 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 572 ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 573 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 574 ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 575 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 576 mult->skipNumeric = PETSC_TRUE; /* a numeric product is done here */ 577 mult->destroy = AB->ops->destroy; 578 mult->duplicate = AB->ops->duplicate; 579 AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32; 580 AB->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_32; 581 AB->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_32; 582 AB->ops->matmult = MatMatMult_MPIAIJ_MPIAIJ; 583 *C = AB; 584 PetscFunctionReturn(0); 585 } 586 587 #undef __FUNCT__ 588 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 589 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 590 { 591 PetscErrorCode ierr; 592 593 PetscFunctionBegin; 594 if (scall == MAT_INITIAL_MATRIX){ 595 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 596 } 597 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 598 PetscFunctionReturn(0); 599 } 600 601 typedef struct { 602 Mat workB; 603 PetscScalar *rvalues,*svalues; 604 MPI_Request *rwaits,*swaits; 605 } MPIAIJ_MPIDense; 606 607 #undef __FUNCT__ 608 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 609 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 610 { 611 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 612 PetscErrorCode ierr; 613 614 PetscFunctionBegin; 615 ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 616 ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 617 ierr = PetscFree(contents);CHKERRQ(ierr); 618 PetscFunctionReturn(0); 619 } 620 621 #undef __FUNCT__ 622 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 623 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 624 { 625 PetscErrorCode ierr; 626 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 627 PetscInt nz = aij->B->cmap->n; 628 PetscContainer container; 629 MPIAIJ_MPIDense *contents; 630 VecScatter ctx = aij->Mvctx; 631 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 632 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 633 PetscInt m=A->rmap->n,n=B->cmap->n; 634 635 PetscFunctionBegin; 636 ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 637 ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 638 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 639 ierr = MatMPIDenseSetPreallocation(*C,PETSC_NULL);CHKERRQ(ierr); 640 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 641 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 642 (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 643 644 ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 645 /* Create work matrix used to store off processor rows of B needed for local product */ 646 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 647 /* Create work arrays needed */ 648 ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 649 B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 650 from->n,MPI_Request,&contents->rwaits, 651 to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 652 653 ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 654 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 655 ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 656 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 657 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 658 PetscFunctionReturn(0); 659 } 660 661 #undef __FUNCT__ 662 #define __FUNCT__ "MatMPIDenseScatter" 663 /* 664 Performs an efficient scatter on the rows of B needed by this process; this is 665 a modification of the VecScatterBegin_() routines. 666 */ 667 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 668 { 669 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 670 PetscErrorCode ierr; 671 PetscScalar *b,*w,*svalues,*rvalues; 672 VecScatter ctx = aij->Mvctx; 673 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 674 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 675 PetscInt i,j,k; 676 PetscInt *sindices,*sstarts,*rindices,*rstarts; 677 PetscMPIInt *sprocs,*rprocs,nrecvs; 678 MPI_Request *swaits,*rwaits; 679 MPI_Comm comm = ((PetscObject)A)->comm; 680 PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 681 MPI_Status status; 682 MPIAIJ_MPIDense *contents; 683 PetscContainer container; 684 Mat workB; 685 686 PetscFunctionBegin; 687 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 688 if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 689 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 690 691 workB = *outworkB = contents->workB; 692 if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n); 693 sindices = to->indices; 694 sstarts = to->starts; 695 sprocs = to->procs; 696 swaits = contents->swaits; 697 svalues = contents->svalues; 698 699 rindices = from->indices; 700 rstarts = from->starts; 701 rprocs = from->procs; 702 rwaits = contents->rwaits; 703 rvalues = contents->rvalues; 704 705 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 706 ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 707 708 for (i=0; i<from->n; i++) { 709 ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 710 } 711 712 for (i=0; i<to->n; i++) { 713 /* pack a message at a time */ 714 CHKMEMQ; 715 for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 716 for (k=0; k<ncols; k++) { 717 svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 718 } 719 } 720 CHKMEMQ; 721 ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 722 } 723 724 nrecvs = from->n; 725 while (nrecvs) { 726 ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 727 nrecvs--; 728 /* unpack a message at a time */ 729 CHKMEMQ; 730 for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 731 for (k=0; k<ncols; k++) { 732 w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 733 } 734 } 735 CHKMEMQ; 736 } 737 if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 738 739 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 740 ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 741 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 742 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 743 PetscFunctionReturn(0); 744 } 745 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 749 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 750 { 751 PetscErrorCode ierr; 752 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 753 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 754 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 755 Mat workB; 756 757 PetscFunctionBegin; 758 759 /* diagonal block of A times all local rows of B*/ 760 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 761 762 /* get off processor parts of B needed to complete the product */ 763 ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 764 765 /* off-diagonal block of A times nonlocal rows of B */ 766 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 767 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 768 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 769 PetscFunctionReturn(0); 770 } 771 772 #undef __FUNCT__ 773 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 774 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C) 775 { 776 PetscErrorCode ierr; 777 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 778 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 779 Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 780 PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 781 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 782 Mat_SeqAIJ *p_loc,*p_oth; 783 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 784 PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 785 PetscInt cm=C->rmap->n,anz,pnz; 786 Mat_PtAPMPI *ptap=c->ptap; 787 PetscScalar *apa_sparse=ptap->apa; 788 PetscInt *api,*apj,*apJ,i,j,k,row; 789 PetscInt cstart=C->cmap->rstart; 790 PetscInt cdnz,conz,k0,k1,nextp; 791 #if defined(DEBUG_MATMATMULT) 792 PetscMPIInt rank=a->rank; 793 #endif 794 795 PetscFunctionBegin; 796 #if defined(DEBUG_MATMATMULT) 797 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 798 #endif 799 800 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 801 /*-----------------------------------------------------*/ 802 /* update numerical values of P_oth and P_loc */ 803 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 804 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 805 #if defined(DEBUG_MATMATMULT) 806 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 807 #endif 808 809 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 810 /*----------------------------------------------------------*/ 811 /* get data from symbolic products */ 812 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 813 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 814 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 815 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 816 817 api = ptap->api; 818 apj = ptap->apj; 819 for (i=0; i<cm; i++) { 820 apJ = apj + api[i]; 821 822 /* diagonal portion of A */ 823 anz = adi[i+1] - adi[i]; 824 adj = ad->j + adi[i]; 825 ada = ad->a + adi[i]; 826 for (j=0; j<anz; j++) { 827 row = adj[j]; 828 pnz = pi_loc[row+1] - pi_loc[row]; 829 pj = pj_loc + pi_loc[row]; 830 pa = pa_loc + pi_loc[row]; 831 /* perform sparse axpy */ 832 valtmp = ada[j]; 833 nextp = 0; 834 for (k=0; nextp<pnz; k++) { 835 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 836 apa_sparse[k] += valtmp*pa[nextp++]; 837 } 838 } 839 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 840 } 841 842 /* off-diagonal portion of A */ 843 anz = aoi[i+1] - aoi[i]; 844 aoj = ao->j + aoi[i]; 845 aoa = ao->a + aoi[i]; 846 for (j=0; j<anz; j++) { 847 row = aoj[j]; 848 pnz = pi_oth[row+1] - pi_oth[row]; 849 pj = pj_oth + pi_oth[row]; 850 pa = pa_oth + pi_oth[row]; 851 /* perform sparse axpy */ 852 valtmp = aoa[j]; 853 nextp = 0; 854 for (k=0; nextp<pnz; k++) { 855 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 856 apa_sparse[k] += valtmp*pa[nextp++]; 857 } 858 } 859 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 860 } 861 862 /* set values in C */ 863 cdnz = cd->i[i+1] - cd->i[i]; 864 conz = co->i[i+1] - co->i[i]; 865 866 /* 1st off-diagoanl part of C */ 867 ca = coa + co->i[i]; 868 k = 0; 869 for (k0=0; k0<conz; k0++){ 870 if (apJ[k] >= cstart) break; 871 ca[k0] = apa_sparse[k]; 872 apa_sparse[k] = 0.0; 873 k++; 874 } 875 876 /* diagonal part of C */ 877 ca = cda + cd->i[i]; 878 for (k1=0; k1<cdnz; k1++){ 879 ca[k1] = apa_sparse[k]; 880 apa_sparse[k] = 0.0; 881 k++; 882 } 883 884 /* 2nd off-diagoanl part of C */ 885 ca = coa + co->i[i]; 886 for (; k0<conz; k0++){ 887 ca[k0] = apa_sparse[k]; 888 apa_sparse[k] = 0.0; 889 k++; 890 } 891 } 892 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 893 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 894 PetscFunctionReturn(0); 895 } 896 897 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */ 898 #undef __FUNCT__ 899 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 900 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C) 901 { 902 PetscErrorCode ierr; 903 MPI_Comm comm=((PetscObject)A)->comm; 904 Mat Cmpi; 905 Mat_PtAPMPI *ptap; 906 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 907 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 908 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 909 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 910 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 911 PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 912 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 913 PetscInt nlnk_max,armax,prmax; 914 PetscReal afill; 915 PetscScalar *apa; 916 #if defined(DEBUG_MATMATMULT) 917 PetscMPIInt rank=a->rank; 918 #endif 919 920 PetscFunctionBegin; 921 #if defined(DEBUG_MATMATMULT) 922 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 923 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 924 #endif 925 926 /* create struct Mat_PtAPMPI and attached it to C later */ 927 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 928 929 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 930 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 931 #if defined(DEBUG_MATMATMULT) 932 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 933 #endif 934 /* get P_loc by taking all local rows of P */ 935 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 936 937 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 938 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 939 pi_loc = p_loc->i; pj_loc = p_loc->j; 940 pi_oth = p_oth->i; pj_oth = p_oth->j; 941 942 #if defined(DEBUG_MATMATMULT) 943 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); 944 #endif 945 946 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 947 /*-------------------------------------------------------------------*/ 948 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 949 ptap->api = api; 950 api[0] = 0; 951 952 /* create and initialize a linked list */ 953 armax = ad->rmax+ao->rmax; 954 prmax = PetscMax(p_loc->rmax,p_oth->rmax); 955 nlnk_max = armax*prmax; 956 if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 957 #if defined(DEBUG_MATMATMULT) 958 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] pN %d; nlnk_max %d; Armax %d+%d=%d; Prmax Max(%d,%d)=%d\n",rank,pN,nlnk_max,ad->rmax,ao->rmax,armax,p_loc->rmax,p_oth->rmax,prmax); 959 #endif 960 ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 961 962 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 963 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 964 current_space = free_space; 965 966 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 967 for (i=0; i<am; i++) { 968 apnz = 0; 969 /* diagonal portion of A */ 970 nzi = adi[i+1] - adi[i]; 971 for (j=0; j<nzi; j++){ 972 row = *adj++; 973 pnz = pi_loc[row+1] - pi_loc[row]; 974 Jptr = pj_loc + pi_loc[row]; 975 /* add non-zero cols of P into the sorted linked list lnk */ 976 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 977 } 978 /* off-diagonal portion of A */ 979 nzi = aoi[i+1] - aoi[i]; 980 for (j=0; j<nzi; j++){ 981 row = *aoj++; 982 pnz = pi_oth[row+1] - pi_oth[row]; 983 Jptr = pj_oth + pi_oth[row]; 984 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 985 } 986 987 apnz = *lnk; 988 api[i+1] = api[i] + apnz; 989 if (apnz > apnz_max) apnz_max = apnz; 990 991 /* if free space is not available, double the total space in the list */ 992 if (current_space->local_remaining<apnz) { 993 ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 994 nspacedouble++; 995 } 996 997 /* Copy data into free space, then initialize lnk */ 998 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 999 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1000 current_space->array += apnz; 1001 current_space->local_used += apnz; 1002 current_space->local_remaining -= apnz; 1003 } 1004 1005 /* Allocate space for apj, initialize apj, and */ 1006 /* destroy list of free space and other temporary array(s) */ 1007 ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 1008 apj = ptap->apj; 1009 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 1010 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 1011 #if defined(DEBUG_MATMATMULT) 1012 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done..., apnz_max %d\n",rank,apnz_max); 1013 #endif 1014 1015 /* create and assemble symbolic parallel matrix Cmpi */ 1016 /*----------------------------------------------------*/ 1017 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1018 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1019 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1020 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1021 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1022 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1023 1024 /* malloc apa for assembly Cmpi */ 1025 ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 1026 ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr); 1027 ptap->apa = apa; 1028 for (i=0; i<am; i++){ 1029 row = i + rstart; 1030 apnz = api[i+1] - api[i]; 1031 ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 1032 apj += apnz; 1033 } 1034 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1035 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1036 1037 ptap->destroy = Cmpi->ops->destroy; 1038 ptap->duplicate = Cmpi->ops->duplicate; 1039 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 1040 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1041 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 1042 1043 /* attach the supporting struct to Cmpi for reuse */ 1044 c = (Mat_MPIAIJ*)Cmpi->data; 1045 c->ptap = ptap; 1046 1047 *C = Cmpi; 1048 1049 /* set MatInfo */ 1050 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 1051 if (afill < 1.0) afill = 1.0; 1052 Cmpi->info.mallocs = nspacedouble; 1053 Cmpi->info.fill_ratio_given = fill; 1054 Cmpi->info.fill_ratio_needed = afill; 1055 1056 #if defined(PETSC_USE_INFO) 1057 if (api[am]) { 1058 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 1059 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 1060 } else { 1061 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1062 } 1063 #endif 1064 PetscFunctionReturn(0); 1065 } 1066 1067 /*-------------------------------------------------------------------------*/ 1068 #undef __FUNCT__ 1069 #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ" 1070 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 1071 { 1072 PetscErrorCode ierr; 1073 PetscBool scalable=PETSC_FALSE; 1074 1075 PetscFunctionBegin; 1076 if (scall == MAT_INITIAL_MATRIX){ 1077 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 1078 ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 1079 if (scalable){ 1080 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr); 1081 } else { 1082 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1083 } 1084 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1085 } 1086 ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 1087 PetscFunctionReturn(0); 1088 } 1089 1090 #undef __FUNCT__ 1091 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ" 1092 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1093 { 1094 PetscErrorCode ierr; 1095 Mat_Merge_SeqsToMPI *merge; 1096 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1097 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1098 Mat_PtAPMPI *ptap; 1099 PetscInt *adj,*aJ; 1100 PetscInt i,j,k,anz,pnz,row,*cj; 1101 MatScalar *ada,*aval,*ca,valtmp; 1102 PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1103 MPI_Comm comm=((PetscObject)C)->comm; 1104 PetscMPIInt size,rank,taga,*len_s; 1105 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1106 PetscInt **buf_ri,**buf_rj; 1107 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1108 MPI_Request *s_waits,*r_waits; 1109 MPI_Status *status; 1110 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1111 PetscInt *ai,*aj,*coi,*coj; 1112 PetscInt *poJ=po->j,*pdJ=pd->j; 1113 Mat A_loc; 1114 Mat_SeqAIJ *a_loc; 1115 1116 PetscFunctionBegin; 1117 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1118 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1119 1120 ptap = c->ptap; 1121 merge = ptap->merge; 1122 1123 #if defined(DEBUG_MATTrMatMult) 1124 ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultNumeric: Crmax %d \n",rank,ptap->rmax); 1125 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1126 #endif 1127 1128 /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1129 /*--------------------------------------------------------------*/ 1130 /* get data from symbolic products */ 1131 coi = merge->coi; coj = merge->coj; 1132 ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 1133 ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 1134 1135 bi = merge->bi; bj = merge->bj; 1136 owners = merge->rowmap->range; 1137 ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1138 ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1139 1140 /* get A_loc by taking all local rows of A */ 1141 A_loc = ptap->A_loc; 1142 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1143 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1144 ai = a_loc->i; 1145 aj = a_loc->j; 1146 1147 ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */ 1148 ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 1149 1150 for (i=0; i<am; i++) { 1151 /* 2-a) put A[i,:] to dense array aval */ 1152 anz = ai[i+1] - ai[i]; 1153 adj = aj + ai[i]; 1154 ada = a_loc->a + ai[i]; 1155 for (j=0; j<anz; j++){ 1156 aval[adj[j]] = ada[j]; 1157 } 1158 1159 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1160 /*--------------------------------------------------------------*/ 1161 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1162 pnz = po->i[i+1] - po->i[i]; 1163 poJ = po->j + po->i[i]; 1164 pA = po->a + po->i[i]; 1165 for (j=0; j<pnz; j++){ 1166 row = poJ[j]; 1167 cnz = coi[row+1] - coi[row]; 1168 cj = coj + coi[row]; 1169 ca = coa + coi[row]; 1170 /* perform dense axpy */ 1171 valtmp = pA[j]; 1172 for (k=0; k<cnz; k++) { 1173 ca[k] += valtmp*aval[cj[k]]; 1174 } 1175 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1176 } 1177 1178 /* put the value into Cd (diagonal part) */ 1179 pnz = pd->i[i+1] - pd->i[i]; 1180 pdJ = pd->j + pd->i[i]; 1181 pA = pd->a + pd->i[i]; 1182 for (j=0; j<pnz; j++){ 1183 row = pdJ[j]; 1184 cnz = bi[row+1] - bi[row]; 1185 cj = bj + bi[row]; 1186 ca = ba + bi[row]; 1187 /* perform dense axpy */ 1188 valtmp = pA[j]; 1189 for (k=0; k<cnz; k++) { 1190 ca[k] += valtmp*aval[cj[k]]; 1191 } 1192 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1193 } 1194 1195 /* zero the current row of Pt*A */ 1196 aJ = aj + ai[i]; 1197 for (k=0; k<anz; k++) aval[aJ[k]] = 0.0; 1198 } 1199 1200 /* 3) send and recv matrix values coa */ 1201 /*------------------------------------*/ 1202 buf_ri = merge->buf_ri; 1203 buf_rj = merge->buf_rj; 1204 len_s = merge->len_s; 1205 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1206 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1207 1208 ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr); 1209 for (proc=0,k=0; proc<size; proc++){ 1210 if (!len_s[proc]) continue; 1211 i = merge->owners_co[proc]; 1212 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1213 k++; 1214 } 1215 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1216 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1217 1218 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1219 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1220 ierr = PetscFree(coa);CHKERRQ(ierr); 1221 1222 /* 4) insert local Cseq and received values into Cmpi */ 1223 /*----------------------------------------------------*/ 1224 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1225 for (k=0; k<merge->nrecv; k++){ 1226 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1227 nrows = *(buf_ri_k[k]); 1228 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1229 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1230 } 1231 1232 for (i=0; i<cm; i++) { 1233 row = owners[rank] + i; /* global row index of C_seq */ 1234 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1235 ba_i = ba + bi[i]; 1236 bnz = bi[i+1] - bi[i]; 1237 /* add received vals into ba */ 1238 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1239 /* i-th row */ 1240 if (i == *nextrow[k]) { 1241 cnz = *(nextci[k]+1) - *nextci[k]; 1242 cj = buf_rj[k] + *(nextci[k]); 1243 ca = abuf_r[k] + *(nextci[k]); 1244 nextcj = 0; 1245 for (j=0; nextcj<cnz; j++){ 1246 if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 1247 ba_i[j] += ca[nextcj++]; 1248 } 1249 } 1250 nextrow[k]++; nextci[k]++; 1251 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1252 } 1253 } 1254 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1255 } 1256 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1257 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1258 1259 ierr = PetscFree(ba);CHKERRQ(ierr); 1260 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1261 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1262 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1263 ierr = PetscFree(aval);CHKERRQ(ierr); 1264 PetscFunctionReturn(0); 1265 } 1266 1267 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1268 #undef __FUNCT__ 1269 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ" 1270 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1271 { 1272 PetscErrorCode ierr; 1273 Mat Cmpi,A_loc,POt,PDt; 1274 Mat_PtAPMPI *ptap; 1275 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1276 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1277 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1278 PetscInt nnz; 1279 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1280 PetscInt am=A->rmap->n,pn=P->cmap->n; 1281 PetscBT lnkbt; 1282 MPI_Comm comm=((PetscObject)A)->comm; 1283 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1284 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1285 PetscInt len,proc,*dnz,*onz,*owners; 1286 PetscInt nzi,*bi,*bj; 1287 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1288 MPI_Request *swaits,*rwaits; 1289 MPI_Status *sstatus,rstatus; 1290 Mat_Merge_SeqsToMPI *merge; 1291 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1292 PetscReal afill=1.0,afill_tmp; 1293 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1294 PetscScalar *vals; 1295 Mat_SeqAIJ *a_loc, *pdt,*pot; 1296 1297 PetscFunctionBegin; 1298 /* check if matrix local sizes are compatible */ 1299 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){ 1300 SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 1301 } 1302 1303 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1304 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1305 #if defined(DEBUG_MATTrMatMult) 1306 ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultSymbolic P: %d %d, %d %d; A %d %d, %d %d\n",rank,P->rmap->N,P->cmap->N,P->rmap->n,P->cmap->n,A->rmap->N,aN,A->rmap->n,A->cmap->n); 1307 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1308 #endif 1309 1310 /* create struct Mat_PtAPMPI and attached it to C later */ 1311 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 1312 1313 /* get A_loc by taking all local rows of A */ 1314 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1315 ptap->A_loc = A_loc; 1316 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1317 ai = a_loc->i; 1318 aj = a_loc->j; 1319 1320 /* determine symbolic Co=(p->B)^T*A - send to others */ 1321 /*----------------------------------------------------*/ 1322 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1323 pdt = (Mat_SeqAIJ*)PDt->data; 1324 pdti = pdt->i; pdtj = pdt->j; 1325 1326 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1327 pot = (Mat_SeqAIJ*)POt->data; 1328 poti = pot->i; potj = pot->j; 1329 1330 /* then, compute symbolic Co = (p->B)^T*A */ 1331 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1332 >= (num of nonzero rows of C_seq) - pn */ 1333 ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 1334 coi[0] = 0; 1335 1336 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1337 nnz = fill*(poti[pon] + ai[am]); 1338 ierr = PetscFreeSpaceGet(nnz,&free_space); 1339 current_space = free_space; 1340 #if defined(DEBUG_MATTrMatMult) 1341 ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],ai[am]);CHKERRQ(ierr); 1342 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1343 #endif 1344 /* create and initialize a linked list */ 1345 i = PetscMax(pdt->rmax,pot->rmax); 1346 Crmax = i*a_loc->rmax*size; 1347 if (!Crmax || Crmax > aN) Crmax = aN; 1348 #if defined(DEBUG_MATTrMatMult) 1349 printf("[%d] rmax A_loc %d * max(PD %d, PO %d)=%d, Crmax %d\n",rank,a_loc->rmax,pdt->rmax,pot->rmax,i*a_loc->rmax,Crmax); 1350 #endif 1351 ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1352 1353 for (i=0; i<pon; i++) { 1354 pnz = poti[i+1] - poti[i]; 1355 ptJ = potj + poti[i]; 1356 for (j=0; j<pnz; j++){ 1357 row = ptJ[j]; /* row of A_loc == col of Pot */ 1358 anz = ai[row+1] - ai[row]; 1359 Jptr = aj + ai[row]; 1360 /* add non-zero cols of AP into the sorted linked list lnk */ 1361 ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1362 } 1363 nnz = lnk[0]; 1364 1365 /* If free space is not available, double the total space in the list */ 1366 if (current_space->local_remaining<nnz) { 1367 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1368 nspacedouble++; 1369 } 1370 1371 /* Copy data into free space, and zero out denserows */ 1372 ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1373 current_space->array += nnz; 1374 current_space->local_used += nnz; 1375 current_space->local_remaining -= nnz; 1376 coi[i+1] = coi[i] + nnz; 1377 } 1378 1379 ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 1380 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1381 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]); 1382 if (afill_tmp > afill) afill = afill_tmp; 1383 1384 /* send j-array (coj) of Co to other processors */ 1385 /*----------------------------------------------*/ 1386 /* determine row ownership */ 1387 ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 1388 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1389 merge->rowmap->n = pn; 1390 merge->rowmap->bs = 1; 1391 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1392 owners = merge->rowmap->range; 1393 1394 /* determine the number of messages to send, their lengths */ 1395 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 1396 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1397 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 1398 len_s = merge->len_s; 1399 merge->nsend = 0; 1400 1401 ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 1402 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1403 1404 proc = 0; 1405 for (i=0; i<pon; i++){ 1406 while (prmap[i] >= owners[proc+1]) proc++; 1407 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1408 len_s[proc] += coi[i+1] - coi[i]; 1409 } 1410 1411 len = 0; /* max length of buf_si[] */ 1412 owners_co[0] = 0; 1413 for (proc=0; proc<size; proc++){ 1414 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1415 if (len_si[proc]){ 1416 merge->nsend++; 1417 len_si[proc] = 2*(len_si[proc] + 1); 1418 len += len_si[proc]; 1419 } 1420 } 1421 1422 /* determine the number and length of messages to receive for coi and coj */ 1423 ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1424 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1425 1426 /* post the Irecv and Isend of coj */ 1427 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1428 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1429 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 1430 for (proc=0, k=0; proc<size; proc++){ 1431 if (!len_s[proc]) continue; 1432 i = owners_co[proc]; 1433 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1434 k++; 1435 } 1436 1437 /* receives and sends of coj are complete */ 1438 ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 1439 for (i=0; i<merge->nrecv; i++){ 1440 PetscMPIInt icompleted; 1441 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1442 } 1443 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1444 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1445 1446 /* send and recv coi */ 1447 /*-------------------*/ 1448 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1449 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1450 ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 1451 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1452 for (proc=0,k=0; proc<size; proc++){ 1453 if (!len_s[proc]) continue; 1454 /* form outgoing message for i-structure: 1455 buf_si[0]: nrows to be sent 1456 [1:nrows]: row index (global) 1457 [nrows+1:2*nrows+1]: i-structure index 1458 */ 1459 /*-------------------------------------------*/ 1460 nrows = len_si[proc]/2 - 1; 1461 buf_si_i = buf_si + nrows+1; 1462 buf_si[0] = nrows; 1463 buf_si_i[0] = 0; 1464 nrows = 0; 1465 for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 1466 nzi = coi[i+1] - coi[i]; 1467 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1468 buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 1469 nrows++; 1470 } 1471 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1472 k++; 1473 buf_si += len_si[proc]; 1474 } 1475 i = merge->nrecv; 1476 while (i--) { 1477 PetscMPIInt icompleted; 1478 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1479 } 1480 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1481 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1482 ierr = PetscFree(len_si);CHKERRQ(ierr); 1483 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1484 ierr = PetscFree(swaits);CHKERRQ(ierr); 1485 ierr = PetscFree(sstatus);CHKERRQ(ierr); 1486 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1487 1488 /* compute the local portion of C (mpi mat) */ 1489 /*------------------------------------------*/ 1490 /* allocate bi array and free space for accumulating nonzero column info */ 1491 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1492 bi[0] = 0; 1493 1494 /* set initial free space to be fill*(nnz(P) + nnz(A)) */ 1495 nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1496 ierr = PetscFreeSpaceGet(nnz,&free_space); 1497 current_space = free_space; 1498 1499 #if defined(DEBUG_MATTrMatMult) 1500 ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=fill %g*(%d + %d + %d)\n",rank,fill,pdti[pn],poti[pon],ai[am]); 1501 #endif 1502 1503 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1504 for (k=0; k<merge->nrecv; k++){ 1505 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1506 nrows = *buf_ri_k[k]; 1507 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1508 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1509 } 1510 1511 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1512 rmax = 0; 1513 for (i=0; i<pn; i++) { 1514 /* add pdt[i,:]*AP into lnk */ 1515 pnz = pdti[i+1] - pdti[i]; 1516 ptJ = pdtj + pdti[i]; 1517 for (j=0; j<pnz; j++){ 1518 row = ptJ[j]; /* row of AP == col of Pt */ 1519 anz = ai[row+1] - ai[row]; 1520 Jptr = aj + ai[row]; 1521 /* add non-zero cols of AP into the sorted linked list lnk */ 1522 ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1523 } 1524 1525 /* add received col data into lnk */ 1526 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1527 if (i == *nextrow[k]) { /* i-th row */ 1528 nzi = *(nextci[k]+1) - *nextci[k]; 1529 Jptr = buf_rj[k] + *nextci[k]; 1530 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1531 nextrow[k]++; nextci[k]++; 1532 } 1533 } 1534 nnz = lnk[0]; 1535 1536 /* if free space is not available, make more free space */ 1537 if (current_space->local_remaining<nnz) { 1538 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1539 nspacedouble++; 1540 } 1541 /* copy data into free space, then initialize lnk */ 1542 ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1543 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1544 current_space->array += nnz; 1545 current_space->local_used += nnz; 1546 current_space->local_remaining -= nnz; 1547 bi[i+1] = bi[i] + nnz; 1548 if (nnz > rmax) rmax = nnz; 1549 } 1550 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1551 1552 ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1553 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1554 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]); 1555 if (afill_tmp > afill) afill = afill_tmp; 1556 ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 1557 ierr = MatDestroy(&POt);CHKERRQ(ierr); 1558 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1559 1560 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 1561 /*----------------------------------------------------------------------------------*/ 1562 ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 1563 ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr); 1564 1565 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1566 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1567 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1568 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1569 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1570 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1571 for (i=0; i<pn; i++){ 1572 row = i + rstart; 1573 nnz = bi[i+1] - bi[i]; 1574 Jptr = bj + bi[i]; 1575 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1576 } 1577 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1578 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1579 ierr = PetscFree(vals);CHKERRQ(ierr); 1580 1581 merge->bi = bi; 1582 merge->bj = bj; 1583 merge->coi = coi; 1584 merge->coj = coj; 1585 merge->buf_ri = buf_ri; 1586 merge->buf_rj = buf_rj; 1587 merge->owners_co = owners_co; 1588 merge->destroy = Cmpi->ops->destroy; 1589 merge->duplicate = Cmpi->ops->duplicate; 1590 1591 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 1592 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1593 1594 /* attach the supporting struct to Cmpi for reuse */ 1595 c = (Mat_MPIAIJ*)Cmpi->data; 1596 c->ptap = ptap; 1597 ptap->api = PETSC_NULL; 1598 ptap->apj = PETSC_NULL; 1599 ptap->merge = merge; 1600 ptap->rmax = rmax; 1601 1602 *C = Cmpi; 1603 #if defined(PETSC_USE_INFO) 1604 if (bi[pn] != 0) { 1605 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 1606 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 1607 } else { 1608 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1609 } 1610 #endif 1611 PetscFunctionReturn(0); 1612 } 1613 1614 #undef __FUNCT__ 1615 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 1616 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C) 1617 { 1618 PetscErrorCode ierr; 1619 Mat_Merge_SeqsToMPI *merge; 1620 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1621 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1622 Mat_PtAPMPI *ptap; 1623 PetscInt *adj; 1624 PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1625 MatScalar *ada,*ca,valtmp; 1626 PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1627 MPI_Comm comm=((PetscObject)C)->comm; 1628 PetscMPIInt size,rank,taga,*len_s; 1629 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1630 PetscInt **buf_ri,**buf_rj; 1631 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1632 MPI_Request *s_waits,*r_waits; 1633 MPI_Status *status; 1634 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1635 PetscInt *ai,*aj,*coi,*coj; 1636 PetscInt *poJ=po->j,*pdJ=pd->j; 1637 Mat A_loc; 1638 Mat_SeqAIJ *a_loc; 1639 1640 PetscFunctionBegin; 1641 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1642 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1643 1644 ptap = c->ptap; 1645 merge = ptap->merge; 1646 1647 #if defined(DEBUG_MATTrMatMult) 1648 ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultNumeric_Scalable: Crmax %d \n",rank,ptap->rmax); 1649 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1650 #endif 1651 1652 /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1653 /*------------------------------------------*/ 1654 /* get data from symbolic products */ 1655 coi = merge->coi; coj = merge->coj; 1656 ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 1657 ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 1658 1659 bi = merge->bi; bj = merge->bj; 1660 owners = merge->rowmap->range; 1661 ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1662 ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1663 1664 /* get A_loc by taking all local rows of A */ 1665 A_loc = ptap->A_loc; 1666 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1667 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1668 ai = a_loc->i; 1669 aj = a_loc->j; 1670 1671 for (i=0; i<am; i++) { 1672 /* 2-a) put A[i,:] to dense array aval */ 1673 anz = ai[i+1] - ai[i]; 1674 adj = aj + ai[i]; 1675 ada = a_loc->a + ai[i]; 1676 1677 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1678 /*-------------------------------------------------------------*/ 1679 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1680 pnz = po->i[i+1] - po->i[i]; 1681 poJ = po->j + po->i[i]; 1682 pA = po->a + po->i[i]; 1683 for (j=0; j<pnz; j++){ 1684 row = poJ[j]; 1685 cnz = coi[row+1] - coi[row]; 1686 cj = coj + coi[row]; 1687 ca = coa + coi[row]; 1688 /* perform sparse axpy */ 1689 nexta = 0; 1690 valtmp = pA[j]; 1691 for (k=0; nexta<anz; k++) { 1692 if (cj[k] == adj[nexta]){ 1693 ca[k] += valtmp*ada[nexta]; 1694 nexta++; 1695 } 1696 } 1697 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1698 } 1699 1700 /* put the value into Cd (diagonal part) */ 1701 pnz = pd->i[i+1] - pd->i[i]; 1702 pdJ = pd->j + pd->i[i]; 1703 pA = pd->a + pd->i[i]; 1704 for (j=0; j<pnz; j++){ 1705 row = pdJ[j]; 1706 cnz = bi[row+1] - bi[row]; 1707 cj = bj + bi[row]; 1708 ca = ba + bi[row]; 1709 /* perform sparse axpy */ 1710 nexta = 0; 1711 valtmp = pA[j]; 1712 for (k=0; nexta<anz; k++) { 1713 if (cj[k] == adj[nexta]){ 1714 ca[k] += valtmp*ada[nexta]; 1715 nexta++; 1716 } 1717 } 1718 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1719 } 1720 1721 } 1722 1723 /* 3) send and recv matrix values coa */ 1724 /*------------------------------------*/ 1725 buf_ri = merge->buf_ri; 1726 buf_rj = merge->buf_rj; 1727 len_s = merge->len_s; 1728 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1729 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1730 1731 ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr); 1732 for (proc=0,k=0; proc<size; proc++){ 1733 if (!len_s[proc]) continue; 1734 i = merge->owners_co[proc]; 1735 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1736 k++; 1737 } 1738 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1739 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1740 1741 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1742 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1743 ierr = PetscFree(coa);CHKERRQ(ierr); 1744 1745 /* 4) insert local Cseq and received values into Cmpi */ 1746 /*----------------------------------------------------*/ 1747 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1748 for (k=0; k<merge->nrecv; k++){ 1749 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1750 nrows = *(buf_ri_k[k]); 1751 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1752 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1753 } 1754 1755 for (i=0; i<cm; i++) { 1756 row = owners[rank] + i; /* global row index of C_seq */ 1757 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1758 ba_i = ba + bi[i]; 1759 bnz = bi[i+1] - bi[i]; 1760 /* add received vals into ba */ 1761 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1762 /* i-th row */ 1763 if (i == *nextrow[k]) { 1764 cnz = *(nextci[k]+1) - *nextci[k]; 1765 cj = buf_rj[k] + *(nextci[k]); 1766 ca = abuf_r[k] + *(nextci[k]); 1767 nextcj = 0; 1768 for (j=0; nextcj<cnz; j++){ 1769 if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 1770 ba_i[j] += ca[nextcj++]; 1771 } 1772 } 1773 nextrow[k]++; nextci[k]++; 1774 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1775 } 1776 } 1777 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1778 } 1779 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1780 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1781 1782 ierr = PetscFree(ba);CHKERRQ(ierr); 1783 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1784 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1785 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1786 PetscFunctionReturn(0); 1787 } 1788 1789 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1790 #undef __FUNCT__ 1791 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 1792 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C) 1793 { 1794 PetscErrorCode ierr; 1795 Mat Cmpi,A_loc,POt,PDt; 1796 Mat_PtAPMPI *ptap; 1797 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1798 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1799 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1800 PetscInt nnz; 1801 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1802 PetscInt am=A->rmap->n,pn=P->cmap->n; 1803 MPI_Comm comm=((PetscObject)A)->comm; 1804 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1805 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1806 PetscInt len,proc,*dnz,*onz,*owners; 1807 PetscInt nzi,*bi,*bj; 1808 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1809 MPI_Request *swaits,*rwaits; 1810 MPI_Status *sstatus,rstatus; 1811 Mat_Merge_SeqsToMPI *merge; 1812 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1813 PetscReal afill=1.0,afill_tmp; 1814 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1815 PetscScalar *vals; 1816 Mat_SeqAIJ *a_loc, *pdt,*pot; 1817 1818 PetscFunctionBegin; 1819 /* check if matrix local sizes are compatible */ 1820 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){ 1821 SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 1822 } 1823 1824 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1825 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1826 #if defined(DEBUG_MATTrMatMult) 1827 ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultSymbolic_Scalable P: %d %d, %d %d; A %d %d, %d %d\n",rank,P->rmap->N,P->cmap->N,P->rmap->n,P->cmap->n,A->rmap->N,aN,A->rmap->n,A->cmap->n); 1828 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1829 #endif 1830 1831 /* create struct Mat_PtAPMPI and attached it to C later */ 1832 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 1833 1834 /* get A_loc by taking all local rows of A */ 1835 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1836 ptap->A_loc = A_loc; 1837 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1838 ai = a_loc->i; 1839 aj = a_loc->j; 1840 1841 /* determine symbolic Co=(p->B)^T*A - send to others */ 1842 /*----------------------------------------------------*/ 1843 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1844 pdt = (Mat_SeqAIJ*)PDt->data; 1845 pdti = pdt->i; pdtj = pdt->j; 1846 1847 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1848 pot = (Mat_SeqAIJ*)POt->data; 1849 poti = pot->i; potj = pot->j; 1850 1851 /* then, compute symbolic Co = (p->B)^T*A */ 1852 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1853 >= (num of nonzero rows of C_seq) - pn */ 1854 ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 1855 coi[0] = 0; 1856 1857 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1858 nnz = fill*(poti[pon] + ai[am]); 1859 ierr = PetscFreeSpaceGet(nnz,&free_space); 1860 current_space = free_space; 1861 #if defined(DEBUG_MATTrMatMult) 1862 ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],ai[am]);CHKERRQ(ierr); 1863 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1864 #endif 1865 /* create and initialize a linked list */ 1866 i = PetscMax(pdt->rmax,pot->rmax); 1867 Crmax = i*a_loc->rmax*size; /* non-scalable! */ 1868 if (!Crmax || Crmax > aN) Crmax = aN; 1869 #if defined(DEBUG_MATTrMatMult) 1870 printf("[%d] rmax A_loc %d * max(PD %d, PO %d)=%d, Crmax %d\n",rank,a_loc->rmax,pdt->rmax,pot->rmax,i*a_loc->rmax,Crmax); 1871 #endif 1872 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 1873 1874 for (i=0; i<pon; i++) { 1875 nnz = 0; 1876 pnz = poti[i+1] - poti[i]; 1877 ptJ = potj + poti[i]; 1878 for (j=0; j<pnz; j++){ 1879 row = ptJ[j]; /* row of A_loc == col of Pot */ 1880 anz = ai[row+1] - ai[row]; 1881 Jptr = aj + ai[row]; 1882 /* add non-zero cols of AP into the sorted linked list lnk */ 1883 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1884 } 1885 nnz = lnk[0]; 1886 1887 /* If free space is not available, double the total space in the list */ 1888 if (current_space->local_remaining<nnz) { 1889 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1890 nspacedouble++; 1891 } 1892 1893 /* Copy data into free space, and zero out denserows */ 1894 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1895 current_space->array += nnz; 1896 current_space->local_used += nnz; 1897 current_space->local_remaining -= nnz; 1898 coi[i+1] = coi[i] + nnz; 1899 } 1900 1901 ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 1902 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1903 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]); 1904 if (afill_tmp > afill) afill = afill_tmp; 1905 1906 /* send j-array (coj) of Co to other processors */ 1907 /*----------------------------------------------*/ 1908 /* determine row ownership */ 1909 ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 1910 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1911 merge->rowmap->n = pn; 1912 merge->rowmap->bs = 1; 1913 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1914 owners = merge->rowmap->range; 1915 1916 /* determine the number of messages to send, their lengths */ 1917 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 1918 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1919 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 1920 len_s = merge->len_s; 1921 merge->nsend = 0; 1922 1923 ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 1924 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1925 1926 proc = 0; 1927 for (i=0; i<pon; i++){ 1928 while (prmap[i] >= owners[proc+1]) proc++; 1929 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1930 len_s[proc] += coi[i+1] - coi[i]; 1931 } 1932 1933 len = 0; /* max length of buf_si[] */ 1934 owners_co[0] = 0; 1935 for (proc=0; proc<size; proc++){ 1936 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1937 if (len_si[proc]){ 1938 merge->nsend++; 1939 len_si[proc] = 2*(len_si[proc] + 1); 1940 len += len_si[proc]; 1941 } 1942 } 1943 1944 /* determine the number and length of messages to receive for coi and coj */ 1945 ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1946 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1947 1948 /* post the Irecv and Isend of coj */ 1949 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1950 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1951 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 1952 for (proc=0, k=0; proc<size; proc++){ 1953 if (!len_s[proc]) continue; 1954 i = owners_co[proc]; 1955 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1956 k++; 1957 } 1958 1959 /* receives and sends of coj are complete */ 1960 ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 1961 for (i=0; i<merge->nrecv; i++){ 1962 PetscMPIInt icompleted; 1963 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1964 } 1965 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1966 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1967 1968 /* send and recv coi */ 1969 /*-------------------*/ 1970 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1971 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1972 ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 1973 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1974 for (proc=0,k=0; proc<size; proc++){ 1975 if (!len_s[proc]) continue; 1976 /* form outgoing message for i-structure: 1977 buf_si[0]: nrows to be sent 1978 [1:nrows]: row index (global) 1979 [nrows+1:2*nrows+1]: i-structure index 1980 */ 1981 /*-------------------------------------------*/ 1982 nrows = len_si[proc]/2 - 1; 1983 buf_si_i = buf_si + nrows+1; 1984 buf_si[0] = nrows; 1985 buf_si_i[0] = 0; 1986 nrows = 0; 1987 for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 1988 nzi = coi[i+1] - coi[i]; 1989 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1990 buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 1991 nrows++; 1992 } 1993 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1994 k++; 1995 buf_si += len_si[proc]; 1996 } 1997 i = merge->nrecv; 1998 while (i--) { 1999 PetscMPIInt icompleted; 2000 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 2001 } 2002 ierr = PetscFree(rwaits);CHKERRQ(ierr); 2003 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 2004 ierr = PetscFree(len_si);CHKERRQ(ierr); 2005 ierr = PetscFree(len_ri);CHKERRQ(ierr); 2006 ierr = PetscFree(swaits);CHKERRQ(ierr); 2007 ierr = PetscFree(sstatus);CHKERRQ(ierr); 2008 ierr = PetscFree(buf_s);CHKERRQ(ierr); 2009 2010 /* compute the local portion of C (mpi mat) */ 2011 /*------------------------------------------*/ 2012 /* allocate bi array and free space for accumulating nonzero column info */ 2013 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 2014 bi[0] = 0; 2015 2016 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 2017 nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 2018 ierr = PetscFreeSpaceGet(nnz,&free_space); 2019 current_space = free_space; 2020 2021 #if defined(DEBUG_MATTrMatMult) 2022 ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=%d + %d + %d\n",rank,pdti[pn],poti[pon],ai[am]); 2023 #endif 2024 2025 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 2026 for (k=0; k<merge->nrecv; k++){ 2027 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 2028 nrows = *buf_ri_k[k]; 2029 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 2030 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 2031 } 2032 2033 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 2034 rmax = 0; 2035 for (i=0; i<pn; i++) { 2036 /* add pdt[i,:]*AP into lnk */ 2037 pnz = pdti[i+1] - pdti[i]; 2038 ptJ = pdtj + pdti[i]; 2039 for (j=0; j<pnz; j++){ 2040 row = ptJ[j]; /* row of AP == col of Pt */ 2041 anz = ai[row+1] - ai[row]; 2042 Jptr = aj + ai[row]; 2043 /* add non-zero cols of AP into the sorted linked list lnk */ 2044 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 2045 } 2046 2047 /* add received col data into lnk */ 2048 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 2049 if (i == *nextrow[k]) { /* i-th row */ 2050 nzi = *(nextci[k]+1) - *nextci[k]; 2051 Jptr = buf_rj[k] + *nextci[k]; 2052 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 2053 nextrow[k]++; nextci[k]++; 2054 } 2055 } 2056 nnz = lnk[0]; 2057 2058 /* if free space is not available, make more free space */ 2059 if (current_space->local_remaining<nnz) { 2060 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 2061 nspacedouble++; 2062 } 2063 /* copy data into free space, then initialize lnk */ 2064 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2065 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2066 current_space->array += nnz; 2067 current_space->local_used += nnz; 2068 current_space->local_remaining -= nnz; 2069 bi[i+1] = bi[i] + nnz; 2070 if (nnz > rmax) rmax = nnz; 2071 } 2072 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 2073 2074 ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 2075 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2076 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]); 2077 if (afill_tmp > afill) afill = afill_tmp; 2078 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 2079 ierr = MatDestroy(&POt);CHKERRQ(ierr); 2080 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 2081 2082 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 2083 /*----------------------------------------------------------------------------------*/ 2084 ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 2085 ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr); 2086 2087 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 2088 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 2089 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 2090 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 2091 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2092 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 2093 for (i=0; i<pn; i++){ 2094 row = i + rstart; 2095 nnz = bi[i+1] - bi[i]; 2096 Jptr = bj + bi[i]; 2097 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 2098 } 2099 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2100 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2101 ierr = PetscFree(vals);CHKERRQ(ierr); 2102 2103 merge->bi = bi; 2104 merge->bj = bj; 2105 merge->coi = coi; 2106 merge->coj = coj; 2107 merge->buf_ri = buf_ri; 2108 merge->buf_rj = buf_rj; 2109 merge->owners_co = owners_co; 2110 merge->destroy = Cmpi->ops->destroy; 2111 merge->duplicate = Cmpi->ops->duplicate; 2112 2113 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 2114 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 2115 2116 /* attach the supporting struct to Cmpi for reuse */ 2117 c = (Mat_MPIAIJ*)Cmpi->data; 2118 c->ptap = ptap; 2119 ptap->api = PETSC_NULL; 2120 ptap->apj = PETSC_NULL; 2121 ptap->merge = merge; 2122 ptap->rmax = rmax; 2123 ptap->apa = PETSC_NULL; 2124 2125 *C = Cmpi; 2126 #if defined(PETSC_USE_INFO) 2127 if (bi[pn] != 0) { 2128 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 2129 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 2130 } else { 2131 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 2132 } 2133 #endif 2134 PetscFunctionReturn(0); 2135 } 2136