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