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