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 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_32; 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 PetscFunctionReturn(0); 584 } 585 586 #undef __FUNCT__ 587 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 588 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 589 { 590 PetscErrorCode ierr; 591 592 PetscFunctionBegin; 593 if (scall == MAT_INITIAL_MATRIX){ 594 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 595 } 596 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 597 PetscFunctionReturn(0); 598 } 599 600 typedef struct { 601 Mat workB; 602 PetscScalar *rvalues,*svalues; 603 MPI_Request *rwaits,*swaits; 604 } MPIAIJ_MPIDense; 605 606 #undef __FUNCT__ 607 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 608 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 609 { 610 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 611 PetscErrorCode ierr; 612 613 PetscFunctionBegin; 614 ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 615 ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 616 ierr = PetscFree(contents);CHKERRQ(ierr); 617 PetscFunctionReturn(0); 618 } 619 620 #undef __FUNCT__ 621 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 622 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 623 { 624 PetscErrorCode ierr; 625 Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 626 PetscInt nz = aij->B->cmap->n; 627 PetscContainer container; 628 MPIAIJ_MPIDense *contents; 629 VecScatter ctx = aij->Mvctx; 630 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 631 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 632 PetscInt m=A->rmap->n,n=B->cmap->n; 633 634 PetscFunctionBegin; 635 ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 636 ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 637 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 638 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 639 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 640 (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 641 642 ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 643 /* Create work matrix used to store off processor rows of B needed for local product */ 644 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 645 /* Create work arrays needed */ 646 ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 647 B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 648 from->n,MPI_Request,&contents->rwaits, 649 to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 650 651 ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 652 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 653 ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 654 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 655 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 656 PetscFunctionReturn(0); 657 } 658 659 #undef __FUNCT__ 660 #define __FUNCT__ "MatMPIDenseScatter" 661 /* 662 Performs an efficient scatter on the rows of B needed by this process; this is 663 a modification of the VecScatterBegin_() routines. 664 */ 665 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 666 { 667 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 668 PetscErrorCode ierr; 669 PetscScalar *b,*w,*svalues,*rvalues; 670 VecScatter ctx = aij->Mvctx; 671 VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 672 VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 673 PetscInt i,j,k; 674 PetscInt *sindices,*sstarts,*rindices,*rstarts; 675 PetscMPIInt *sprocs,*rprocs,nrecvs; 676 MPI_Request *swaits,*rwaits; 677 MPI_Comm comm = ((PetscObject)A)->comm; 678 PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 679 MPI_Status status; 680 MPIAIJ_MPIDense *contents; 681 PetscContainer container; 682 Mat workB; 683 684 PetscFunctionBegin; 685 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 686 if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 687 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 688 689 workB = *outworkB = contents->workB; 690 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); 691 sindices = to->indices; 692 sstarts = to->starts; 693 sprocs = to->procs; 694 swaits = contents->swaits; 695 svalues = contents->svalues; 696 697 rindices = from->indices; 698 rstarts = from->starts; 699 rprocs = from->procs; 700 rwaits = contents->rwaits; 701 rvalues = contents->rvalues; 702 703 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 704 ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 705 706 for (i=0; i<from->n; i++) { 707 ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 708 } 709 710 for (i=0; i<to->n; i++) { 711 /* pack a message at a time */ 712 CHKMEMQ; 713 for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 714 for (k=0; k<ncols; k++) { 715 svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 716 } 717 } 718 CHKMEMQ; 719 ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 720 } 721 722 nrecvs = from->n; 723 while (nrecvs) { 724 ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 725 nrecvs--; 726 /* unpack a message at a time */ 727 CHKMEMQ; 728 for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 729 for (k=0; k<ncols; k++) { 730 w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 731 } 732 } 733 CHKMEMQ; 734 } 735 if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 736 737 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 738 ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 739 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 740 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 741 PetscFunctionReturn(0); 742 } 743 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 744 745 #undef __FUNCT__ 746 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 747 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 748 { 749 PetscErrorCode ierr; 750 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 751 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 752 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 753 Mat workB; 754 755 PetscFunctionBegin; 756 757 /* diagonal block of A times all local rows of B*/ 758 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 759 760 /* get off processor parts of B needed to complete the product */ 761 ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 762 763 /* off-diagonal block of A times nonlocal rows of B */ 764 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 765 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 766 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 767 PetscFunctionReturn(0); 768 } 769 770 #undef __FUNCT__ 771 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 772 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C) 773 { 774 PetscErrorCode ierr; 775 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 776 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 777 Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 778 PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 779 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 780 Mat_SeqAIJ *p_loc,*p_oth; 781 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 782 PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 783 PetscInt cm=C->rmap->n,anz,pnz; 784 Mat_PtAPMPI *ptap=c->ptap; 785 PetscScalar *apa_sparse=ptap->apa; 786 PetscInt *api,*apj,*apJ,i,j,k,row; 787 PetscInt cstart=C->cmap->rstart; 788 PetscInt cdnz,conz,k0,k1,nextp; 789 #if defined(DEBUG_MATMATMULT) 790 PetscMPIInt rank=a->rank; 791 #endif 792 793 PetscFunctionBegin; 794 #if defined(DEBUG_MATMATMULT) 795 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 796 #endif 797 798 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 799 /*-----------------------------------------------------*/ 800 /* update numerical values of P_oth and P_loc */ 801 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 802 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 803 #if defined(DEBUG_MATMATMULT) 804 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 805 #endif 806 807 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 808 /*----------------------------------------------------------*/ 809 /* get data from symbolic products */ 810 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 811 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 812 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 813 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 814 815 api = ptap->api; 816 apj = ptap->apj; 817 for (i=0; i<cm; i++) { 818 apJ = apj + api[i]; 819 820 /* diagonal portion of A */ 821 anz = adi[i+1] - adi[i]; 822 adj = ad->j + adi[i]; 823 ada = ad->a + adi[i]; 824 for (j=0; j<anz; j++) { 825 row = adj[j]; 826 pnz = pi_loc[row+1] - pi_loc[row]; 827 pj = pj_loc + pi_loc[row]; 828 pa = pa_loc + pi_loc[row]; 829 /* perform sparse axpy */ 830 valtmp = ada[j]; 831 nextp = 0; 832 for (k=0; nextp<pnz; k++) { 833 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 834 apa_sparse[k] += valtmp*pa[nextp++]; 835 } 836 } 837 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 838 } 839 840 /* off-diagonal portion of A */ 841 anz = aoi[i+1] - aoi[i]; 842 aoj = ao->j + aoi[i]; 843 aoa = ao->a + aoi[i]; 844 for (j=0; j<anz; j++) { 845 row = aoj[j]; 846 pnz = pi_oth[row+1] - pi_oth[row]; 847 pj = pj_oth + pi_oth[row]; 848 pa = pa_oth + pi_oth[row]; 849 /* perform sparse axpy */ 850 valtmp = aoa[j]; 851 nextp = 0; 852 for (k=0; nextp<pnz; k++) { 853 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 854 apa_sparse[k] += valtmp*pa[nextp++]; 855 } 856 } 857 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 858 } 859 860 /* set values in C */ 861 cdnz = cd->i[i+1] - cd->i[i]; 862 conz = co->i[i+1] - co->i[i]; 863 864 /* 1st off-diagoanl part of C */ 865 ca = coa + co->i[i]; 866 k = 0; 867 for (k0=0; k0<conz; k0++){ 868 if (apJ[k] >= cstart) break; 869 ca[k0] = apa_sparse[k]; 870 apa_sparse[k] = 0.0; 871 k++; 872 } 873 874 /* diagonal part of C */ 875 ca = cda + cd->i[i]; 876 for (k1=0; k1<cdnz; k1++){ 877 ca[k1] = apa_sparse[k]; 878 apa_sparse[k] = 0.0; 879 k++; 880 } 881 882 /* 2nd off-diagoanl part of C */ 883 ca = coa + co->i[i]; 884 for (; k0<conz; k0++){ 885 ca[k0] = apa_sparse[k]; 886 apa_sparse[k] = 0.0; 887 k++; 888 } 889 } 890 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 891 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 892 PetscFunctionReturn(0); 893 } 894 895 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */ 896 #undef __FUNCT__ 897 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 898 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C) 899 { 900 PetscErrorCode ierr; 901 MPI_Comm comm=((PetscObject)A)->comm; 902 Mat Cmpi; 903 Mat_PtAPMPI *ptap; 904 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 905 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 906 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 907 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 908 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 909 PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 910 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 911 PetscInt nlnk_max,armax,prmax; 912 PetscReal afill; 913 PetscScalar *apa; 914 #if defined(DEBUG_MATMATMULT) 915 PetscMPIInt rank=a->rank; 916 #endif 917 918 PetscFunctionBegin; 919 #if defined(DEBUG_MATMATMULT) 920 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 921 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 922 #endif 923 924 /* create struct Mat_PtAPMPI and attached it to C later */ 925 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 926 927 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 928 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 929 #if defined(DEBUG_MATMATMULT) 930 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 931 #endif 932 /* get P_loc by taking all local rows of P */ 933 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 934 935 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 936 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 937 pi_loc = p_loc->i; pj_loc = p_loc->j; 938 pi_oth = p_oth->i; pj_oth = p_oth->j; 939 940 #if defined(DEBUG_MATMATMULT) 941 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); 942 #endif 943 944 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 945 /*-------------------------------------------------------------------*/ 946 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 947 ptap->api = api; 948 api[0] = 0; 949 950 /* create and initialize a linked list */ 951 armax = ad->rmax+ao->rmax; 952 prmax = PetscMax(p_loc->rmax,p_oth->rmax); 953 nlnk_max = armax*prmax; 954 if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 955 #if defined(DEBUG_MATMATMULT) 956 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); 957 #endif 958 ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 959 960 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 961 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 962 current_space = free_space; 963 964 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 965 for (i=0; i<am; i++) { 966 apnz = 0; 967 /* diagonal portion of A */ 968 nzi = adi[i+1] - adi[i]; 969 for (j=0; j<nzi; j++){ 970 row = *adj++; 971 pnz = pi_loc[row+1] - pi_loc[row]; 972 Jptr = pj_loc + pi_loc[row]; 973 /* add non-zero cols of P into the sorted linked list lnk */ 974 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 975 } 976 /* off-diagonal portion of A */ 977 nzi = aoi[i+1] - aoi[i]; 978 for (j=0; j<nzi; j++){ 979 row = *aoj++; 980 pnz = pi_oth[row+1] - pi_oth[row]; 981 Jptr = pj_oth + pi_oth[row]; 982 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 983 } 984 985 apnz = *lnk; 986 api[i+1] = api[i] + apnz; 987 if (apnz > apnz_max) apnz_max = apnz; 988 989 /* if free space is not available, double the total space in the list */ 990 if (current_space->local_remaining<apnz) { 991 ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 992 nspacedouble++; 993 } 994 995 /* Copy data into free space, then initialize lnk */ 996 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 997 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 998 current_space->array += apnz; 999 current_space->local_used += apnz; 1000 current_space->local_remaining -= apnz; 1001 } 1002 1003 /* Allocate space for apj, initialize apj, and */ 1004 /* destroy list of free space and other temporary array(s) */ 1005 ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 1006 apj = ptap->apj; 1007 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 1008 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 1009 #if defined(DEBUG_MATMATMULT) 1010 if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done..., apnz_max %d\n",rank,apnz_max); 1011 #endif 1012 1013 /* create and assemble symbolic parallel matrix Cmpi */ 1014 /*----------------------------------------------------*/ 1015 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1016 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1017 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1018 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1019 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1020 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1021 1022 /* malloc apa for assembly Cmpi */ 1023 ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 1024 ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr); 1025 ptap->apa = apa; 1026 for (i=0; i<am; i++){ 1027 row = i + rstart; 1028 apnz = api[i+1] - api[i]; 1029 ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 1030 apj += apnz; 1031 } 1032 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1033 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1034 1035 ptap->destroy = Cmpi->ops->destroy; 1036 ptap->duplicate = Cmpi->ops->duplicate; 1037 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 1038 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1039 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 1040 1041 /* attach the supporting struct to Cmpi for reuse */ 1042 c = (Mat_MPIAIJ*)Cmpi->data; 1043 c->ptap = ptap; 1044 1045 *C = Cmpi; 1046 1047 /* set MatInfo */ 1048 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 1049 if (afill < 1.0) afill = 1.0; 1050 Cmpi->info.mallocs = nspacedouble; 1051 Cmpi->info.fill_ratio_given = fill; 1052 Cmpi->info.fill_ratio_needed = afill; 1053 1054 #if defined(PETSC_USE_INFO) 1055 if (api[am]) { 1056 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 1057 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 1058 } else { 1059 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1060 } 1061 #endif 1062 PetscFunctionReturn(0); 1063 } 1064 1065 /*-------------------------------------------------------------------------*/ 1066 #undef __FUNCT__ 1067 #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ" 1068 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 1069 { 1070 PetscErrorCode ierr; 1071 PetscBool scalable=PETSC_FALSE; 1072 1073 PetscFunctionBegin; 1074 if (scall == MAT_INITIAL_MATRIX){ 1075 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 1076 ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 1077 if (scalable){ 1078 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr); 1079 } else { 1080 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1081 } 1082 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1083 } 1084 ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 1085 PetscFunctionReturn(0); 1086 } 1087 1088 #undef __FUNCT__ 1089 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ" 1090 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1091 { 1092 PetscErrorCode ierr; 1093 Mat_Merge_SeqsToMPI *merge; 1094 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1095 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1096 Mat_PtAPMPI *ptap; 1097 PetscInt *adj,*aJ; 1098 PetscInt i,j,k,anz,pnz,row,*cj; 1099 MatScalar *ada,*aval,*ca,valtmp; 1100 PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1101 MPI_Comm comm=((PetscObject)C)->comm; 1102 PetscMPIInt size,rank,taga,*len_s; 1103 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1104 PetscInt **buf_ri,**buf_rj; 1105 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1106 MPI_Request *s_waits,*r_waits; 1107 MPI_Status *status; 1108 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1109 PetscInt *ai,*aj,*coi,*coj; 1110 PetscInt *poJ=po->j,*pdJ=pd->j; 1111 Mat A_loc; 1112 Mat_SeqAIJ *a_loc; 1113 1114 PetscFunctionBegin; 1115 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1116 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1117 1118 ptap = c->ptap; 1119 merge = ptap->merge; 1120 1121 #if defined(DEBUG_MATTrMatMult) 1122 ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultNumeric: Crmax %d \n",rank,ptap->rmax); 1123 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1124 #endif 1125 1126 /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1127 /*--------------------------------------------------------------*/ 1128 /* get data from symbolic products */ 1129 coi = merge->coi; coj = merge->coj; 1130 ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 1131 ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 1132 1133 bi = merge->bi; bj = merge->bj; 1134 owners = merge->rowmap->range; 1135 ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1136 ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1137 1138 /* get A_loc by taking all local rows of A */ 1139 A_loc = ptap->A_loc; 1140 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1141 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1142 ai = a_loc->i; 1143 aj = a_loc->j; 1144 1145 ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */ 1146 ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 1147 1148 for (i=0; i<am; i++) { 1149 /* 2-a) put A[i,:] to dense array aval */ 1150 anz = ai[i+1] - ai[i]; 1151 adj = aj + ai[i]; 1152 ada = a_loc->a + ai[i]; 1153 for (j=0; j<anz; j++){ 1154 aval[adj[j]] = ada[j]; 1155 } 1156 1157 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1158 /*--------------------------------------------------------------*/ 1159 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1160 pnz = po->i[i+1] - po->i[i]; 1161 poJ = po->j + po->i[i]; 1162 pA = po->a + po->i[i]; 1163 for (j=0; j<pnz; j++){ 1164 row = poJ[j]; 1165 cnz = coi[row+1] - coi[row]; 1166 cj = coj + coi[row]; 1167 ca = coa + coi[row]; 1168 /* perform dense axpy */ 1169 valtmp = pA[j]; 1170 for (k=0; k<cnz; k++) { 1171 ca[k] += valtmp*aval[cj[k]]; 1172 } 1173 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1174 } 1175 1176 /* put the value into Cd (diagonal part) */ 1177 pnz = pd->i[i+1] - pd->i[i]; 1178 pdJ = pd->j + pd->i[i]; 1179 pA = pd->a + pd->i[i]; 1180 for (j=0; j<pnz; j++){ 1181 row = pdJ[j]; 1182 cnz = bi[row+1] - bi[row]; 1183 cj = bj + bi[row]; 1184 ca = ba + bi[row]; 1185 /* perform dense axpy */ 1186 valtmp = pA[j]; 1187 for (k=0; k<cnz; k++) { 1188 ca[k] += valtmp*aval[cj[k]]; 1189 } 1190 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1191 } 1192 1193 /* zero the current row of Pt*A */ 1194 aJ = aj + ai[i]; 1195 for (k=0; k<anz; k++) aval[aJ[k]] = 0.0; 1196 } 1197 1198 /* 3) send and recv matrix values coa */ 1199 /*------------------------------------*/ 1200 buf_ri = merge->buf_ri; 1201 buf_rj = merge->buf_rj; 1202 len_s = merge->len_s; 1203 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1204 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1205 1206 ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr); 1207 for (proc=0,k=0; proc<size; proc++){ 1208 if (!len_s[proc]) continue; 1209 i = merge->owners_co[proc]; 1210 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1211 k++; 1212 } 1213 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1214 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1215 1216 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1217 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1218 ierr = PetscFree(coa);CHKERRQ(ierr); 1219 1220 /* 4) insert local Cseq and received values into Cmpi */ 1221 /*----------------------------------------------------*/ 1222 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1223 for (k=0; k<merge->nrecv; k++){ 1224 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1225 nrows = *(buf_ri_k[k]); 1226 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1227 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1228 } 1229 1230 for (i=0; i<cm; i++) { 1231 row = owners[rank] + i; /* global row index of C_seq */ 1232 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1233 ba_i = ba + bi[i]; 1234 bnz = bi[i+1] - bi[i]; 1235 /* add received vals into ba */ 1236 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1237 /* i-th row */ 1238 if (i == *nextrow[k]) { 1239 cnz = *(nextci[k]+1) - *nextci[k]; 1240 cj = buf_rj[k] + *(nextci[k]); 1241 ca = abuf_r[k] + *(nextci[k]); 1242 nextcj = 0; 1243 for (j=0; nextcj<cnz; j++){ 1244 if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 1245 ba_i[j] += ca[nextcj++]; 1246 } 1247 } 1248 nextrow[k]++; nextci[k]++; 1249 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1250 } 1251 } 1252 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1253 } 1254 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1255 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1256 1257 ierr = PetscFree(ba);CHKERRQ(ierr); 1258 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1259 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1260 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1261 ierr = PetscFree(aval);CHKERRQ(ierr); 1262 PetscFunctionReturn(0); 1263 } 1264 1265 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1266 #undef __FUNCT__ 1267 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ" 1268 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1269 { 1270 PetscErrorCode ierr; 1271 Mat Cmpi,A_loc,POt,PDt; 1272 Mat_PtAPMPI *ptap; 1273 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1274 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1275 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1276 PetscInt nnz; 1277 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1278 PetscInt am=A->rmap->n,pn=P->cmap->n; 1279 PetscBT lnkbt; 1280 MPI_Comm comm=((PetscObject)A)->comm; 1281 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1282 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1283 PetscInt len,proc,*dnz,*onz,*owners; 1284 PetscInt nzi,*bi,*bj; 1285 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1286 MPI_Request *swaits,*rwaits; 1287 MPI_Status *sstatus,rstatus; 1288 Mat_Merge_SeqsToMPI *merge; 1289 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1290 PetscReal afill=1.0,afill_tmp; 1291 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1292 PetscScalar *vals; 1293 Mat_SeqAIJ *a_loc, *pdt,*pot; 1294 1295 PetscFunctionBegin; 1296 /* check if matrix local sizes are compatible */ 1297 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){ 1298 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); 1299 } 1300 1301 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1302 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1303 #if defined(DEBUG_MATTrMatMult) 1304 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); 1305 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1306 #endif 1307 1308 /* create struct Mat_PtAPMPI and attached it to C later */ 1309 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 1310 1311 /* get A_loc by taking all local rows of A */ 1312 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1313 ptap->A_loc = A_loc; 1314 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1315 ai = a_loc->i; 1316 aj = a_loc->j; 1317 1318 /* determine symbolic Co=(p->B)^T*A - send to others */ 1319 /*----------------------------------------------------*/ 1320 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1321 pdt = (Mat_SeqAIJ*)PDt->data; 1322 pdti = pdt->i; pdtj = pdt->j; 1323 1324 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1325 pot = (Mat_SeqAIJ*)POt->data; 1326 poti = pot->i; potj = pot->j; 1327 1328 /* then, compute symbolic Co = (p->B)^T*A */ 1329 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1330 >= (num of nonzero rows of C_seq) - pn */ 1331 ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 1332 coi[0] = 0; 1333 1334 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1335 nnz = fill*(poti[pon] + ai[am]); 1336 ierr = PetscFreeSpaceGet(nnz,&free_space); 1337 current_space = free_space; 1338 #if defined(DEBUG_MATTrMatMult) 1339 ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],ai[am]);CHKERRQ(ierr); 1340 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1341 #endif 1342 /* create and initialize a linked list */ 1343 i = PetscMax(pdt->rmax,pot->rmax); 1344 Crmax = i*a_loc->rmax*size; 1345 if (!Crmax || Crmax > aN) Crmax = aN; 1346 #if defined(DEBUG_MATTrMatMult) 1347 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); 1348 #endif 1349 ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1350 1351 for (i=0; i<pon; i++) { 1352 pnz = poti[i+1] - poti[i]; 1353 ptJ = potj + poti[i]; 1354 for (j=0; j<pnz; j++){ 1355 row = ptJ[j]; /* row of A_loc == col of Pot */ 1356 anz = ai[row+1] - ai[row]; 1357 Jptr = aj + ai[row]; 1358 /* add non-zero cols of AP into the sorted linked list lnk */ 1359 ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1360 } 1361 nnz = lnk[0]; 1362 1363 /* If free space is not available, double the total space in the list */ 1364 if (current_space->local_remaining<nnz) { 1365 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1366 nspacedouble++; 1367 } 1368 1369 /* Copy data into free space, and zero out denserows */ 1370 ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1371 current_space->array += nnz; 1372 current_space->local_used += nnz; 1373 current_space->local_remaining -= nnz; 1374 coi[i+1] = coi[i] + nnz; 1375 } 1376 1377 ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 1378 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1379 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]); 1380 if (afill_tmp > afill) afill = afill_tmp; 1381 1382 /* send j-array (coj) of Co to other processors */ 1383 /*----------------------------------------------*/ 1384 /* determine row ownership */ 1385 ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 1386 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1387 merge->rowmap->n = pn; 1388 merge->rowmap->bs = 1; 1389 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1390 owners = merge->rowmap->range; 1391 1392 /* determine the number of messages to send, their lengths */ 1393 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 1394 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1395 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 1396 len_s = merge->len_s; 1397 merge->nsend = 0; 1398 1399 ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 1400 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1401 1402 proc = 0; 1403 for (i=0; i<pon; i++){ 1404 while (prmap[i] >= owners[proc+1]) proc++; 1405 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1406 len_s[proc] += coi[i+1] - coi[i]; 1407 } 1408 1409 len = 0; /* max length of buf_si[] */ 1410 owners_co[0] = 0; 1411 for (proc=0; proc<size; proc++){ 1412 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1413 if (len_si[proc]){ 1414 merge->nsend++; 1415 len_si[proc] = 2*(len_si[proc] + 1); 1416 len += len_si[proc]; 1417 } 1418 } 1419 1420 /* determine the number and length of messages to receive for coi and coj */ 1421 ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1422 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1423 1424 /* post the Irecv and Isend of coj */ 1425 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1426 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1427 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 1428 for (proc=0, k=0; proc<size; proc++){ 1429 if (!len_s[proc]) continue; 1430 i = owners_co[proc]; 1431 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1432 k++; 1433 } 1434 1435 /* receives and sends of coj are complete */ 1436 ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 1437 for (i=0; i<merge->nrecv; i++){ 1438 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1439 } 1440 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1441 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1442 1443 /* send and recv coi */ 1444 /*-------------------*/ 1445 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1446 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1447 ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 1448 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1449 for (proc=0,k=0; proc<size; proc++){ 1450 if (!len_s[proc]) continue; 1451 /* form outgoing message for i-structure: 1452 buf_si[0]: nrows to be sent 1453 [1:nrows]: row index (global) 1454 [nrows+1:2*nrows+1]: i-structure index 1455 */ 1456 /*-------------------------------------------*/ 1457 nrows = len_si[proc]/2 - 1; 1458 buf_si_i = buf_si + nrows+1; 1459 buf_si[0] = nrows; 1460 buf_si_i[0] = 0; 1461 nrows = 0; 1462 for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 1463 nzi = coi[i+1] - coi[i]; 1464 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1465 buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 1466 nrows++; 1467 } 1468 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1469 k++; 1470 buf_si += len_si[proc]; 1471 } 1472 i = merge->nrecv; 1473 while (i--) { 1474 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1475 } 1476 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1477 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1478 ierr = PetscFree(len_si);CHKERRQ(ierr); 1479 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1480 ierr = PetscFree(swaits);CHKERRQ(ierr); 1481 ierr = PetscFree(sstatus);CHKERRQ(ierr); 1482 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1483 1484 /* compute the local portion of C (mpi mat) */ 1485 /*------------------------------------------*/ 1486 /* allocate bi array and free space for accumulating nonzero column info */ 1487 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1488 bi[0] = 0; 1489 1490 /* set initial free space to be fill*(nnz(P) + nnz(A)) */ 1491 nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1492 ierr = PetscFreeSpaceGet(nnz,&free_space); 1493 current_space = free_space; 1494 1495 #if defined(DEBUG_MATTrMatMult) 1496 ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=fill %g*(%d + %d + %d)\n",rank,fill,pdti[pn],poti[pon],ai[am]); 1497 #endif 1498 1499 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1500 for (k=0; k<merge->nrecv; k++){ 1501 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1502 nrows = *buf_ri_k[k]; 1503 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1504 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1505 } 1506 1507 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1508 rmax = 0; 1509 for (i=0; i<pn; i++) { 1510 /* add pdt[i,:]*AP into lnk */ 1511 pnz = pdti[i+1] - pdti[i]; 1512 ptJ = pdtj + pdti[i]; 1513 for (j=0; j<pnz; j++){ 1514 row = ptJ[j]; /* row of AP == col of Pt */ 1515 anz = ai[row+1] - ai[row]; 1516 Jptr = aj + ai[row]; 1517 /* add non-zero cols of AP into the sorted linked list lnk */ 1518 ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1519 } 1520 1521 /* add received col data into lnk */ 1522 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1523 if (i == *nextrow[k]) { /* i-th row */ 1524 nzi = *(nextci[k]+1) - *nextci[k]; 1525 Jptr = buf_rj[k] + *nextci[k]; 1526 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1527 nextrow[k]++; nextci[k]++; 1528 } 1529 } 1530 nnz = lnk[0]; 1531 1532 /* if free space is not available, make more free space */ 1533 if (current_space->local_remaining<nnz) { 1534 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1535 nspacedouble++; 1536 } 1537 /* copy data into free space, then initialize lnk */ 1538 ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1539 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1540 current_space->array += nnz; 1541 current_space->local_used += nnz; 1542 current_space->local_remaining -= nnz; 1543 bi[i+1] = bi[i] + nnz; 1544 if (nnz > rmax) rmax = nnz; 1545 } 1546 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1547 1548 ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1549 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1550 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]); 1551 if (afill_tmp > afill) afill = afill_tmp; 1552 ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 1553 ierr = MatDestroy(&POt);CHKERRQ(ierr); 1554 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1555 1556 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 1557 /*----------------------------------------------------------------------------------*/ 1558 ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 1559 ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr); 1560 1561 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1562 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1563 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1564 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1565 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1566 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1567 for (i=0; i<pn; i++){ 1568 row = i + rstart; 1569 nnz = bi[i+1] - bi[i]; 1570 Jptr = bj + bi[i]; 1571 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1572 } 1573 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1574 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1575 ierr = PetscFree(vals);CHKERRQ(ierr); 1576 1577 merge->bi = bi; 1578 merge->bj = bj; 1579 merge->coi = coi; 1580 merge->coj = coj; 1581 merge->buf_ri = buf_ri; 1582 merge->buf_rj = buf_rj; 1583 merge->owners_co = owners_co; 1584 merge->destroy = Cmpi->ops->destroy; 1585 merge->duplicate = Cmpi->ops->duplicate; 1586 1587 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 1588 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1589 1590 /* attach the supporting struct to Cmpi for reuse */ 1591 c = (Mat_MPIAIJ*)Cmpi->data; 1592 c->ptap = ptap; 1593 ptap->api = PETSC_NULL; 1594 ptap->apj = PETSC_NULL; 1595 ptap->merge = merge; 1596 ptap->rmax = rmax; 1597 1598 *C = Cmpi; 1599 #if defined(PETSC_USE_INFO) 1600 if (bi[pn] != 0) { 1601 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 1602 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 1603 } else { 1604 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1605 } 1606 #endif 1607 PetscFunctionReturn(0); 1608 } 1609 1610 #undef __FUNCT__ 1611 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 1612 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C) 1613 { 1614 PetscErrorCode ierr; 1615 Mat_Merge_SeqsToMPI *merge; 1616 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1617 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1618 Mat_PtAPMPI *ptap; 1619 PetscInt *adj; 1620 PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1621 MatScalar *ada,*ca,valtmp; 1622 PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1623 MPI_Comm comm=((PetscObject)C)->comm; 1624 PetscMPIInt size,rank,taga,*len_s; 1625 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1626 PetscInt **buf_ri,**buf_rj; 1627 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1628 MPI_Request *s_waits,*r_waits; 1629 MPI_Status *status; 1630 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1631 PetscInt *ai,*aj,*coi,*coj; 1632 PetscInt *poJ=po->j,*pdJ=pd->j; 1633 Mat A_loc; 1634 Mat_SeqAIJ *a_loc; 1635 1636 PetscFunctionBegin; 1637 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1638 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1639 1640 ptap = c->ptap; 1641 merge = ptap->merge; 1642 1643 #if defined(DEBUG_MATTrMatMult) 1644 ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultNumeric_Scalable: Crmax %d \n",rank,ptap->rmax); 1645 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1646 #endif 1647 1648 /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1649 /*------------------------------------------*/ 1650 /* get data from symbolic products */ 1651 coi = merge->coi; coj = merge->coj; 1652 ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 1653 ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 1654 1655 bi = merge->bi; bj = merge->bj; 1656 owners = merge->rowmap->range; 1657 ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1658 ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1659 1660 /* get A_loc by taking all local rows of A */ 1661 A_loc = ptap->A_loc; 1662 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1663 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1664 ai = a_loc->i; 1665 aj = a_loc->j; 1666 1667 for (i=0; i<am; i++) { 1668 /* 2-a) put A[i,:] to dense array aval */ 1669 anz = ai[i+1] - ai[i]; 1670 adj = aj + ai[i]; 1671 ada = a_loc->a + ai[i]; 1672 1673 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1674 /*-------------------------------------------------------------*/ 1675 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1676 pnz = po->i[i+1] - po->i[i]; 1677 poJ = po->j + po->i[i]; 1678 pA = po->a + po->i[i]; 1679 for (j=0; j<pnz; j++){ 1680 row = poJ[j]; 1681 cnz = coi[row+1] - coi[row]; 1682 cj = coj + coi[row]; 1683 ca = coa + coi[row]; 1684 /* perform sparse axpy */ 1685 nexta = 0; 1686 valtmp = pA[j]; 1687 for (k=0; nexta<anz; k++) { 1688 if (cj[k] == adj[nexta]){ 1689 ca[k] += valtmp*ada[nexta]; 1690 nexta++; 1691 } 1692 } 1693 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1694 } 1695 1696 /* put the value into Cd (diagonal part) */ 1697 pnz = pd->i[i+1] - pd->i[i]; 1698 pdJ = pd->j + pd->i[i]; 1699 pA = pd->a + pd->i[i]; 1700 for (j=0; j<pnz; j++){ 1701 row = pdJ[j]; 1702 cnz = bi[row+1] - bi[row]; 1703 cj = bj + bi[row]; 1704 ca = ba + bi[row]; 1705 /* perform sparse axpy */ 1706 nexta = 0; 1707 valtmp = pA[j]; 1708 for (k=0; nexta<anz; k++) { 1709 if (cj[k] == adj[nexta]){ 1710 ca[k] += valtmp*ada[nexta]; 1711 nexta++; 1712 } 1713 } 1714 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1715 } 1716 1717 } 1718 1719 /* 3) send and recv matrix values coa */ 1720 /*------------------------------------*/ 1721 buf_ri = merge->buf_ri; 1722 buf_rj = merge->buf_rj; 1723 len_s = merge->len_s; 1724 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1725 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1726 1727 ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr); 1728 for (proc=0,k=0; proc<size; proc++){ 1729 if (!len_s[proc]) continue; 1730 i = merge->owners_co[proc]; 1731 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1732 k++; 1733 } 1734 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1735 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1736 1737 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1738 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1739 ierr = PetscFree(coa);CHKERRQ(ierr); 1740 1741 /* 4) insert local Cseq and received values into Cmpi */ 1742 /*----------------------------------------------------*/ 1743 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1744 for (k=0; k<merge->nrecv; k++){ 1745 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1746 nrows = *(buf_ri_k[k]); 1747 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1748 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1749 } 1750 1751 for (i=0; i<cm; i++) { 1752 row = owners[rank] + i; /* global row index of C_seq */ 1753 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1754 ba_i = ba + bi[i]; 1755 bnz = bi[i+1] - bi[i]; 1756 /* add received vals into ba */ 1757 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1758 /* i-th row */ 1759 if (i == *nextrow[k]) { 1760 cnz = *(nextci[k]+1) - *nextci[k]; 1761 cj = buf_rj[k] + *(nextci[k]); 1762 ca = abuf_r[k] + *(nextci[k]); 1763 nextcj = 0; 1764 for (j=0; nextcj<cnz; j++){ 1765 if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 1766 ba_i[j] += ca[nextcj++]; 1767 } 1768 } 1769 nextrow[k]++; nextci[k]++; 1770 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1771 } 1772 } 1773 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1774 } 1775 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1776 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1777 1778 ierr = PetscFree(ba);CHKERRQ(ierr); 1779 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1780 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1781 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1782 PetscFunctionReturn(0); 1783 } 1784 1785 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1786 #undef __FUNCT__ 1787 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 1788 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C) 1789 { 1790 PetscErrorCode ierr; 1791 Mat Cmpi,A_loc,POt,PDt; 1792 Mat_PtAPMPI *ptap; 1793 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1794 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1795 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1796 PetscInt nnz; 1797 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1798 PetscInt am=A->rmap->n,pn=P->cmap->n; 1799 MPI_Comm comm=((PetscObject)A)->comm; 1800 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1801 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1802 PetscInt len,proc,*dnz,*onz,*owners; 1803 PetscInt nzi,*bi,*bj; 1804 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1805 MPI_Request *swaits,*rwaits; 1806 MPI_Status *sstatus,rstatus; 1807 Mat_Merge_SeqsToMPI *merge; 1808 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1809 PetscReal afill=1.0,afill_tmp; 1810 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1811 PetscScalar *vals; 1812 Mat_SeqAIJ *a_loc, *pdt,*pot; 1813 1814 PetscFunctionBegin; 1815 /* check if matrix local sizes are compatible */ 1816 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){ 1817 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); 1818 } 1819 1820 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1821 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1822 #if defined(DEBUG_MATTrMatMult) 1823 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); 1824 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1825 #endif 1826 1827 /* create struct Mat_PtAPMPI and attached it to C later */ 1828 ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 1829 1830 /* get A_loc by taking all local rows of A */ 1831 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1832 ptap->A_loc = A_loc; 1833 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1834 ai = a_loc->i; 1835 aj = a_loc->j; 1836 1837 /* determine symbolic Co=(p->B)^T*A - send to others */ 1838 /*----------------------------------------------------*/ 1839 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1840 pdt = (Mat_SeqAIJ*)PDt->data; 1841 pdti = pdt->i; pdtj = pdt->j; 1842 1843 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1844 pot = (Mat_SeqAIJ*)POt->data; 1845 poti = pot->i; potj = pot->j; 1846 1847 /* then, compute symbolic Co = (p->B)^T*A */ 1848 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1849 >= (num of nonzero rows of C_seq) - pn */ 1850 ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 1851 coi[0] = 0; 1852 1853 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1854 nnz = fill*(poti[pon] + ai[am]); 1855 ierr = PetscFreeSpaceGet(nnz,&free_space); 1856 current_space = free_space; 1857 #if defined(DEBUG_MATTrMatMult) 1858 ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],ai[am]);CHKERRQ(ierr); 1859 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1860 #endif 1861 /* create and initialize a linked list */ 1862 i = PetscMax(pdt->rmax,pot->rmax); 1863 Crmax = i*a_loc->rmax*size; /* non-scalable! */ 1864 if (!Crmax || Crmax > aN) Crmax = aN; 1865 #if defined(DEBUG_MATTrMatMult) 1866 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); 1867 #endif 1868 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 1869 1870 for (i=0; i<pon; i++) { 1871 nnz = 0; 1872 pnz = poti[i+1] - poti[i]; 1873 ptJ = potj + poti[i]; 1874 for (j=0; j<pnz; j++){ 1875 row = ptJ[j]; /* row of A_loc == col of Pot */ 1876 anz = ai[row+1] - ai[row]; 1877 Jptr = aj + ai[row]; 1878 /* add non-zero cols of AP into the sorted linked list lnk */ 1879 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1880 } 1881 nnz = lnk[0]; 1882 1883 /* If free space is not available, double the total space in the list */ 1884 if (current_space->local_remaining<nnz) { 1885 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1886 nspacedouble++; 1887 } 1888 1889 /* Copy data into free space, and zero out denserows */ 1890 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1891 current_space->array += nnz; 1892 current_space->local_used += nnz; 1893 current_space->local_remaining -= nnz; 1894 coi[i+1] = coi[i] + nnz; 1895 } 1896 1897 ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 1898 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1899 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]); 1900 if (afill_tmp > afill) afill = afill_tmp; 1901 1902 /* send j-array (coj) of Co to other processors */ 1903 /*----------------------------------------------*/ 1904 /* determine row ownership */ 1905 ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 1906 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1907 merge->rowmap->n = pn; 1908 merge->rowmap->bs = 1; 1909 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1910 owners = merge->rowmap->range; 1911 1912 /* determine the number of messages to send, their lengths */ 1913 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 1914 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1915 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 1916 len_s = merge->len_s; 1917 merge->nsend = 0; 1918 1919 ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 1920 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1921 1922 proc = 0; 1923 for (i=0; i<pon; i++){ 1924 while (prmap[i] >= owners[proc+1]) proc++; 1925 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1926 len_s[proc] += coi[i+1] - coi[i]; 1927 } 1928 1929 len = 0; /* max length of buf_si[] */ 1930 owners_co[0] = 0; 1931 for (proc=0; proc<size; proc++){ 1932 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1933 if (len_si[proc]){ 1934 merge->nsend++; 1935 len_si[proc] = 2*(len_si[proc] + 1); 1936 len += len_si[proc]; 1937 } 1938 } 1939 1940 /* determine the number and length of messages to receive for coi and coj */ 1941 ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1942 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1943 1944 /* post the Irecv and Isend of coj */ 1945 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1946 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1947 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 1948 for (proc=0, k=0; proc<size; proc++){ 1949 if (!len_s[proc]) continue; 1950 i = owners_co[proc]; 1951 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1952 k++; 1953 } 1954 1955 /* receives and sends of coj are complete */ 1956 ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 1957 for (i=0; i<merge->nrecv; i++){ 1958 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1959 } 1960 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1961 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1962 1963 /* send and recv coi */ 1964 /*-------------------*/ 1965 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1966 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1967 ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 1968 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1969 for (proc=0,k=0; proc<size; proc++){ 1970 if (!len_s[proc]) continue; 1971 /* form outgoing message for i-structure: 1972 buf_si[0]: nrows to be sent 1973 [1:nrows]: row index (global) 1974 [nrows+1:2*nrows+1]: i-structure index 1975 */ 1976 /*-------------------------------------------*/ 1977 nrows = len_si[proc]/2 - 1; 1978 buf_si_i = buf_si + nrows+1; 1979 buf_si[0] = nrows; 1980 buf_si_i[0] = 0; 1981 nrows = 0; 1982 for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 1983 nzi = coi[i+1] - coi[i]; 1984 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1985 buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 1986 nrows++; 1987 } 1988 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1989 k++; 1990 buf_si += len_si[proc]; 1991 } 1992 i = merge->nrecv; 1993 while (i--) { 1994 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1995 } 1996 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1997 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1998 ierr = PetscFree(len_si);CHKERRQ(ierr); 1999 ierr = PetscFree(len_ri);CHKERRQ(ierr); 2000 ierr = PetscFree(swaits);CHKERRQ(ierr); 2001 ierr = PetscFree(sstatus);CHKERRQ(ierr); 2002 ierr = PetscFree(buf_s);CHKERRQ(ierr); 2003 2004 /* compute the local portion of C (mpi mat) */ 2005 /*------------------------------------------*/ 2006 /* allocate bi array and free space for accumulating nonzero column info */ 2007 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 2008 bi[0] = 0; 2009 2010 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 2011 nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 2012 ierr = PetscFreeSpaceGet(nnz,&free_space); 2013 current_space = free_space; 2014 2015 #if defined(DEBUG_MATTrMatMult) 2016 ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=%d + %d + %d\n",rank,pdti[pn],poti[pon],ai[am]); 2017 #endif 2018 2019 ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 2020 for (k=0; k<merge->nrecv; k++){ 2021 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 2022 nrows = *buf_ri_k[k]; 2023 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 2024 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 2025 } 2026 2027 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 2028 rmax = 0; 2029 for (i=0; i<pn; i++) { 2030 /* add pdt[i,:]*AP into lnk */ 2031 pnz = pdti[i+1] - pdti[i]; 2032 ptJ = pdtj + pdti[i]; 2033 for (j=0; j<pnz; j++){ 2034 row = ptJ[j]; /* row of AP == col of Pt */ 2035 anz = ai[row+1] - ai[row]; 2036 Jptr = aj + ai[row]; 2037 /* add non-zero cols of AP into the sorted linked list lnk */ 2038 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 2039 } 2040 2041 /* add received col data into lnk */ 2042 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 2043 if (i == *nextrow[k]) { /* i-th row */ 2044 nzi = *(nextci[k]+1) - *nextci[k]; 2045 Jptr = buf_rj[k] + *nextci[k]; 2046 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 2047 nextrow[k]++; nextci[k]++; 2048 } 2049 } 2050 nnz = lnk[0]; 2051 2052 /* if free space is not available, make more free space */ 2053 if (current_space->local_remaining<nnz) { 2054 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 2055 nspacedouble++; 2056 } 2057 /* copy data into free space, then initialize lnk */ 2058 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2059 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2060 current_space->array += nnz; 2061 current_space->local_used += nnz; 2062 current_space->local_remaining -= nnz; 2063 bi[i+1] = bi[i] + nnz; 2064 if (nnz > rmax) rmax = nnz; 2065 } 2066 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 2067 2068 ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 2069 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2070 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]); 2071 if (afill_tmp > afill) afill = afill_tmp; 2072 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 2073 ierr = MatDestroy(&POt);CHKERRQ(ierr); 2074 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 2075 2076 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 2077 /*----------------------------------------------------------------------------------*/ 2078 ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 2079 ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr); 2080 2081 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 2082 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 2083 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 2084 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 2085 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2086 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 2087 for (i=0; i<pn; i++){ 2088 row = i + rstart; 2089 nnz = bi[i+1] - bi[i]; 2090 Jptr = bj + bi[i]; 2091 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 2092 } 2093 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2094 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2095 ierr = PetscFree(vals);CHKERRQ(ierr); 2096 2097 merge->bi = bi; 2098 merge->bj = bj; 2099 merge->coi = coi; 2100 merge->coj = coj; 2101 merge->buf_ri = buf_ri; 2102 merge->buf_rj = buf_rj; 2103 merge->owners_co = owners_co; 2104 merge->destroy = Cmpi->ops->destroy; 2105 merge->duplicate = Cmpi->ops->duplicate; 2106 2107 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 2108 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 2109 2110 /* attach the supporting struct to Cmpi for reuse */ 2111 c = (Mat_MPIAIJ*)Cmpi->data; 2112 c->ptap = ptap; 2113 ptap->api = PETSC_NULL; 2114 ptap->apj = PETSC_NULL; 2115 ptap->merge = merge; 2116 ptap->rmax = rmax; 2117 ptap->apa = PETSC_NULL; 2118 2119 *C = Cmpi; 2120 #if defined(PETSC_USE_INFO) 2121 if (bi[pn] != 0) { 2122 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 2123 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 2124 } else { 2125 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 2126 } 2127 #endif 2128 PetscFunctionReturn(0); 2129 } 2130