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