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