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