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