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