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