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 = p_loc->rmax+p_oth->rmax; 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 PetscInt nlnk_max,armax,prmax; 712 PetscReal afill; 713 PetscScalar *apa; 714 PetscTable ta; 715 716 PetscFunctionBegin; 717 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 718 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 719 720 /* create struct Mat_PtAPMPI and attached it to C later */ 721 ierr = PetscNew(&ptap);CHKERRQ(ierr); 722 723 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 724 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 725 726 /* get P_loc by taking all local rows of P */ 727 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 728 729 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 730 pi_loc = p_loc->i; pj_loc = p_loc->j; 731 if (size > 1) { 732 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 733 pi_oth = p_oth->i; pj_oth = p_oth->j; 734 } else { 735 p_oth = NULL; 736 pi_oth = NULL; pj_oth = NULL; 737 } 738 739 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 740 /*-------------------------------------------------------------------*/ 741 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 742 ptap->api = api; 743 api[0] = 0; 744 745 /* create and initialize a linked list */ 746 Crmax = p_loc->rmax+p_oth->rmax; 747 if (Crmax > pN) Crmax = pN; 748 ierr = PetscTableCreate(Crmax,pN,&ta);CHKERRQ(ierr); 749 MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 750 MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 751 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 752 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 753 754 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 755 756 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 757 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 758 759 current_space = free_space; 760 761 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 762 for (i=0; i<am; i++) { 763 /* diagonal portion of A */ 764 nzi = adi[i+1] - adi[i]; 765 for (j=0; j<nzi; j++) { 766 row = *adj++; 767 pnz = pi_loc[row+1] - pi_loc[row]; 768 Jptr = pj_loc + pi_loc[row]; 769 /* add non-zero cols of P into the sorted linked list lnk */ 770 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 771 } 772 /* off-diagonal portion of A */ 773 nzi = aoi[i+1] - aoi[i]; 774 for (j=0; j<nzi; j++) { 775 row = *aoj++; 776 pnz = pi_oth[row+1] - pi_oth[row]; 777 Jptr = pj_oth + pi_oth[row]; 778 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 779 } 780 781 apnz = *lnk; 782 api[i+1] = api[i] + apnz; 783 if (apnz > apnz_max) apnz_max = apnz; 784 785 /* if free space is not available, double the total space in the list */ 786 if (current_space->local_remaining<apnz) { 787 ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 788 nspacedouble++; 789 } 790 791 /* Copy data into free space, then initialize lnk */ 792 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 793 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 794 795 current_space->array += apnz; 796 current_space->local_used += apnz; 797 current_space->local_remaining -= apnz; 798 } 799 800 /* Allocate space for apj, initialize apj, and */ 801 /* destroy list of free space and other temporary array(s) */ 802 ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 803 apj = ptap->apj; 804 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 805 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 806 807 /* create and assemble symbolic parallel matrix Cmpi */ 808 /*----------------------------------------------------*/ 809 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 810 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 811 ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 812 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 813 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 814 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 815 816 /* malloc apa for assembly Cmpi */ 817 ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr); 818 819 ptap->apa = apa; 820 for (i=0; i<am; i++) { 821 row = i + rstart; 822 apnz = api[i+1] - api[i]; 823 ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 824 apj += apnz; 825 } 826 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 827 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 828 829 ptap->destroy = Cmpi->ops->destroy; 830 ptap->duplicate = Cmpi->ops->duplicate; 831 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 832 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 833 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 834 835 /* attach the supporting struct to Cmpi for reuse */ 836 c = (Mat_MPIAIJ*)Cmpi->data; 837 c->ptap = ptap; 838 839 *C = Cmpi; 840 841 /* set MatInfo */ 842 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 843 if (afill < 1.0) afill = 1.0; 844 Cmpi->info.mallocs = nspacedouble; 845 Cmpi->info.fill_ratio_given = fill; 846 Cmpi->info.fill_ratio_needed = afill; 847 848 #if defined(PETSC_USE_INFO) 849 if (api[am]) { 850 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 851 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 852 } else { 853 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 854 } 855 #endif 856 PetscFunctionReturn(0); 857 } 858 859 /*-------------------------------------------------------------------------*/ 860 #undef __FUNCT__ 861 #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ" 862 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 863 { 864 PetscErrorCode ierr; 865 const char *algTypes[3] = {"scalable","nonscalable","matmatmult"}; 866 PetscInt alg=0; /* set default algorithm */ 867 868 PetscFunctionBegin; 869 if (scall == MAT_INITIAL_MATRIX) { 870 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 871 ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr); 872 ierr = PetscOptionsEnd();CHKERRQ(ierr); 873 874 ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 875 switch (alg) { 876 case 1: 877 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr); 878 break; 879 case 2: 880 { 881 Mat Pt; 882 Mat_PtAPMPI *ptap; 883 Mat_MPIAIJ *c; 884 ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr); 885 ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr); 886 c = (Mat_MPIAIJ*)(*C)->data; 887 ptap = c->ptap; 888 ptap->Pt = Pt; 889 (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult; 890 PetscFunctionReturn(0); 891 } 892 break; 893 default: 894 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 895 break; 896 } 897 ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 898 } 899 ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 900 ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 901 ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 902 PetscFunctionReturn(0); 903 } 904 905 /* This routine only works when scall=MAT_REUSE_MATRIX! */ 906 #undef __FUNCT__ 907 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult" 908 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C) 909 { 910 PetscErrorCode ierr; 911 Mat_MPIAIJ *c=(Mat_MPIAIJ*)C->data; 912 Mat_PtAPMPI *ptap= c->ptap; 913 Mat Pt=ptap->Pt; 914 915 PetscFunctionBegin; 916 ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr); 917 ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr); 918 PetscFunctionReturn(0); 919 } 920 921 /* Non-scalable version, use dense axpy */ 922 #undef __FUNCT__ 923 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable" 924 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C) 925 { 926 PetscErrorCode ierr; 927 Mat_Merge_SeqsToMPI *merge; 928 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 929 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 930 Mat_PtAPMPI *ptap; 931 PetscInt *adj,*aJ; 932 PetscInt i,j,k,anz,pnz,row,*cj; 933 MatScalar *ada,*aval,*ca,valtmp; 934 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 935 MPI_Comm comm; 936 PetscMPIInt size,rank,taga,*len_s; 937 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 938 PetscInt **buf_ri,**buf_rj; 939 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 940 MPI_Request *s_waits,*r_waits; 941 MPI_Status *status; 942 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 943 PetscInt *ai,*aj,*coi,*coj; 944 PetscInt *poJ,*pdJ; 945 Mat A_loc; 946 Mat_SeqAIJ *a_loc; 947 948 PetscFunctionBegin; 949 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 950 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 951 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 952 953 ptap = c->ptap; 954 merge = ptap->merge; 955 956 /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 957 /*--------------------------------------------------------------*/ 958 /* get data from symbolic products */ 959 coi = merge->coi; coj = merge->coj; 960 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 961 962 bi = merge->bi; bj = merge->bj; 963 owners = merge->rowmap->range; 964 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 965 966 /* get A_loc by taking all local rows of A */ 967 A_loc = ptap->A_loc; 968 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 969 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 970 ai = a_loc->i; 971 aj = a_loc->j; 972 973 ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */ 974 975 for (i=0; i<am; i++) { 976 /* 2-a) put A[i,:] to dense array aval */ 977 anz = ai[i+1] - ai[i]; 978 adj = aj + ai[i]; 979 ada = a_loc->a + ai[i]; 980 for (j=0; j<anz; j++) { 981 aval[adj[j]] = ada[j]; 982 } 983 984 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 985 /*--------------------------------------------------------------*/ 986 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 987 pnz = po->i[i+1] - po->i[i]; 988 poJ = po->j + po->i[i]; 989 pA = po->a + po->i[i]; 990 for (j=0; j<pnz; j++) { 991 row = poJ[j]; 992 cnz = coi[row+1] - coi[row]; 993 cj = coj + coi[row]; 994 ca = coa + coi[row]; 995 /* perform dense axpy */ 996 valtmp = pA[j]; 997 for (k=0; k<cnz; k++) { 998 ca[k] += valtmp*aval[cj[k]]; 999 } 1000 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1001 } 1002 1003 /* put the value into Cd (diagonal part) */ 1004 pnz = pd->i[i+1] - pd->i[i]; 1005 pdJ = pd->j + pd->i[i]; 1006 pA = pd->a + pd->i[i]; 1007 for (j=0; j<pnz; j++) { 1008 row = pdJ[j]; 1009 cnz = bi[row+1] - bi[row]; 1010 cj = bj + bi[row]; 1011 ca = ba + bi[row]; 1012 /* perform dense axpy */ 1013 valtmp = pA[j]; 1014 for (k=0; k<cnz; k++) { 1015 ca[k] += valtmp*aval[cj[k]]; 1016 } 1017 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1018 } 1019 1020 /* zero the current row of Pt*A */ 1021 aJ = aj + ai[i]; 1022 for (k=0; k<anz; k++) aval[aJ[k]] = 0.0; 1023 } 1024 1025 /* 3) send and recv matrix values coa */ 1026 /*------------------------------------*/ 1027 buf_ri = merge->buf_ri; 1028 buf_rj = merge->buf_rj; 1029 len_s = merge->len_s; 1030 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1031 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1032 1033 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1034 for (proc=0,k=0; proc<size; proc++) { 1035 if (!len_s[proc]) continue; 1036 i = merge->owners_co[proc]; 1037 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1038 k++; 1039 } 1040 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1041 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1042 1043 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1044 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1045 ierr = PetscFree(coa);CHKERRQ(ierr); 1046 1047 /* 4) insert local Cseq and received values into Cmpi */ 1048 /*----------------------------------------------------*/ 1049 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1050 for (k=0; k<merge->nrecv; k++) { 1051 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1052 nrows = *(buf_ri_k[k]); 1053 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1054 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1055 } 1056 1057 for (i=0; i<cm; i++) { 1058 row = owners[rank] + i; /* global row index of C_seq */ 1059 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1060 ba_i = ba + bi[i]; 1061 bnz = bi[i+1] - bi[i]; 1062 /* add received vals into ba */ 1063 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1064 /* i-th row */ 1065 if (i == *nextrow[k]) { 1066 cnz = *(nextci[k]+1) - *nextci[k]; 1067 cj = buf_rj[k] + *(nextci[k]); 1068 ca = abuf_r[k] + *(nextci[k]); 1069 nextcj = 0; 1070 for (j=0; nextcj<cnz; j++) { 1071 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1072 ba_i[j] += ca[nextcj++]; 1073 } 1074 } 1075 nextrow[k]++; nextci[k]++; 1076 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1077 } 1078 } 1079 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1080 } 1081 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1082 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1083 1084 ierr = PetscFree(ba);CHKERRQ(ierr); 1085 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1086 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1087 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1088 ierr = PetscFree(aval);CHKERRQ(ierr); 1089 PetscFunctionReturn(0); 1090 } 1091 1092 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*); 1093 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1094 #undef __FUNCT__ 1095 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable" 1096 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C) 1097 { 1098 PetscErrorCode ierr; 1099 Mat Cmpi,A_loc,POt,PDt; 1100 Mat_PtAPMPI *ptap; 1101 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1102 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c; 1103 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1104 PetscInt nnz; 1105 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1106 PetscInt am=A->rmap->n,pn=P->cmap->n; 1107 PetscBT lnkbt; 1108 MPI_Comm comm; 1109 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1110 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1111 PetscInt len,proc,*dnz,*onz,*owners; 1112 PetscInt nzi,*bi,*bj; 1113 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1114 MPI_Request *swaits,*rwaits; 1115 MPI_Status *sstatus,rstatus; 1116 Mat_Merge_SeqsToMPI *merge; 1117 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1118 PetscReal afill =1.0,afill_tmp; 1119 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N; 1120 PetscScalar *vals; 1121 Mat_SeqAIJ *a_loc, *pdt,*pot; 1122 1123 PetscFunctionBegin; 1124 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1125 /* check if matrix local sizes are compatible */ 1126 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 1127 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); 1128 } 1129 1130 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1131 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1132 1133 /* create struct Mat_PtAPMPI and attached it to C later */ 1134 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1135 1136 /* get A_loc by taking all local rows of A */ 1137 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1138 1139 ptap->A_loc = A_loc; 1140 1141 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1142 ai = a_loc->i; 1143 aj = a_loc->j; 1144 1145 /* determine symbolic Co=(p->B)^T*A - send to others */ 1146 /*----------------------------------------------------*/ 1147 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1148 pdt = (Mat_SeqAIJ*)PDt->data; 1149 pdti = pdt->i; pdtj = pdt->j; 1150 1151 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1152 pot = (Mat_SeqAIJ*)POt->data; 1153 poti = pot->i; potj = pot->j; 1154 1155 /* then, compute symbolic Co = (p->B)^T*A */ 1156 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */ 1157 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1158 coi[0] = 0; 1159 1160 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1161 nnz = fill*(poti[pon] + ai[am]); 1162 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1163 current_space = free_space; 1164 1165 /* create and initialize a linked list */ 1166 ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1167 1168 for (i=0; i<pon; i++) { 1169 pnz = poti[i+1] - poti[i]; 1170 ptJ = potj + poti[i]; 1171 for (j=0; j<pnz; j++) { 1172 row = ptJ[j]; /* row of A_loc == col of Pot */ 1173 anz = ai[row+1] - ai[row]; 1174 Jptr = aj + ai[row]; 1175 /* add non-zero cols of AP into the sorted linked list lnk */ 1176 ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1177 } 1178 nnz = lnk[0]; 1179 1180 /* If free space is not available, double the total space in the list */ 1181 if (current_space->local_remaining<nnz) { 1182 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1183 nspacedouble++; 1184 } 1185 1186 /* Copy data into free space, and zero out denserows */ 1187 ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1188 1189 current_space->array += nnz; 1190 current_space->local_used += nnz; 1191 current_space->local_remaining -= nnz; 1192 1193 coi[i+1] = coi[i] + nnz; 1194 } 1195 1196 ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1197 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1198 1199 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1200 if (afill_tmp > afill) afill = afill_tmp; 1201 1202 /* send j-array (coj) of Co to other processors */ 1203 /*----------------------------------------------*/ 1204 /* determine row ownership */ 1205 ierr = PetscNew(&merge);CHKERRQ(ierr); 1206 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1207 1208 merge->rowmap->n = pn; 1209 merge->rowmap->bs = 1; 1210 1211 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1212 owners = merge->rowmap->range; 1213 1214 /* determine the number of messages to send, their lengths */ 1215 ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1216 ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 1217 1218 len_s = merge->len_s; 1219 merge->nsend = 0; 1220 1221 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1222 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1223 1224 proc = 0; 1225 for (i=0; i<pon; i++) { 1226 while (prmap[i] >= owners[proc+1]) proc++; 1227 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1228 len_s[proc] += coi[i+1] - coi[i]; 1229 } 1230 1231 len = 0; /* max length of buf_si[] */ 1232 owners_co[0] = 0; 1233 for (proc=0; proc<size; proc++) { 1234 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1235 if (len_si[proc]) { 1236 merge->nsend++; 1237 len_si[proc] = 2*(len_si[proc] + 1); 1238 len += len_si[proc]; 1239 } 1240 } 1241 1242 /* determine the number and length of messages to receive for coi and coj */ 1243 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1244 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1245 1246 /* post the Irecv and Isend of coj */ 1247 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1248 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1249 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1250 for (proc=0, k=0; proc<size; proc++) { 1251 if (!len_s[proc]) continue; 1252 i = owners_co[proc]; 1253 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1254 k++; 1255 } 1256 1257 /* receives and sends of coj are complete */ 1258 ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1259 for (i=0; i<merge->nrecv; i++) { 1260 PetscMPIInt icompleted; 1261 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1262 } 1263 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1264 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1265 1266 /* send and recv coi */ 1267 /*-------------------*/ 1268 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1269 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1270 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1271 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1272 for (proc=0,k=0; proc<size; proc++) { 1273 if (!len_s[proc]) continue; 1274 /* form outgoing message for i-structure: 1275 buf_si[0]: nrows to be sent 1276 [1:nrows]: row index (global) 1277 [nrows+1:2*nrows+1]: i-structure index 1278 */ 1279 /*-------------------------------------------*/ 1280 nrows = len_si[proc]/2 - 1; 1281 buf_si_i = buf_si + nrows+1; 1282 buf_si[0] = nrows; 1283 buf_si_i[0] = 0; 1284 nrows = 0; 1285 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1286 nzi = coi[i+1] - coi[i]; 1287 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1288 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1289 nrows++; 1290 } 1291 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1292 k++; 1293 buf_si += len_si[proc]; 1294 } 1295 i = merge->nrecv; 1296 while (i--) { 1297 PetscMPIInt icompleted; 1298 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1299 } 1300 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1301 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1302 ierr = PetscFree(len_si);CHKERRQ(ierr); 1303 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1304 ierr = PetscFree(swaits);CHKERRQ(ierr); 1305 ierr = PetscFree(sstatus);CHKERRQ(ierr); 1306 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1307 1308 /* compute the local portion of C (mpi mat) */ 1309 /*------------------------------------------*/ 1310 /* allocate bi array and free space for accumulating nonzero column info */ 1311 ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1312 bi[0] = 0; 1313 1314 /* set initial free space to be fill*(nnz(P) + nnz(A)) */ 1315 nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1316 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1317 current_space = free_space; 1318 1319 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1320 for (k=0; k<merge->nrecv; k++) { 1321 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1322 nrows = *buf_ri_k[k]; 1323 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1324 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1325 } 1326 1327 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1328 rmax = 0; 1329 for (i=0; i<pn; i++) { 1330 /* add pdt[i,:]*AP into lnk */ 1331 pnz = pdti[i+1] - pdti[i]; 1332 ptJ = pdtj + pdti[i]; 1333 for (j=0; j<pnz; j++) { 1334 row = ptJ[j]; /* row of AP == col of Pt */ 1335 anz = ai[row+1] - ai[row]; 1336 Jptr = aj + ai[row]; 1337 /* add non-zero cols of AP into the sorted linked list lnk */ 1338 ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1339 } 1340 1341 /* add received col data into lnk */ 1342 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1343 if (i == *nextrow[k]) { /* i-th row */ 1344 nzi = *(nextci[k]+1) - *nextci[k]; 1345 Jptr = buf_rj[k] + *nextci[k]; 1346 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1347 nextrow[k]++; nextci[k]++; 1348 } 1349 } 1350 nnz = lnk[0]; 1351 1352 /* if free space is not available, make more free space */ 1353 if (current_space->local_remaining<nnz) { 1354 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1355 nspacedouble++; 1356 } 1357 /* copy data into free space, then initialize lnk */ 1358 ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1359 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1360 1361 current_space->array += nnz; 1362 current_space->local_used += nnz; 1363 current_space->local_remaining -= nnz; 1364 1365 bi[i+1] = bi[i] + nnz; 1366 if (nnz > rmax) rmax = nnz; 1367 } 1368 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1369 1370 ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 1371 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1372 1373 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 1374 if (afill_tmp > afill) afill = afill_tmp; 1375 ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 1376 ierr = MatDestroy(&POt);CHKERRQ(ierr); 1377 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1378 1379 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 1380 /*----------------------------------------------------------------------------------*/ 1381 ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 1382 1383 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1384 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1385 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 1386 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1387 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1388 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1389 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1390 for (i=0; i<pn; i++) { 1391 row = i + rstart; 1392 nnz = bi[i+1] - bi[i]; 1393 Jptr = bj + bi[i]; 1394 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1395 } 1396 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1397 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1398 ierr = PetscFree(vals);CHKERRQ(ierr); 1399 1400 merge->bi = bi; 1401 merge->bj = bj; 1402 merge->coi = coi; 1403 merge->coj = coj; 1404 merge->buf_ri = buf_ri; 1405 merge->buf_rj = buf_rj; 1406 merge->owners_co = owners_co; 1407 merge->destroy = Cmpi->ops->destroy; 1408 merge->duplicate = Cmpi->ops->duplicate; 1409 1410 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1411 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1412 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1413 1414 /* attach the supporting struct to Cmpi for reuse */ 1415 c = (Mat_MPIAIJ*)Cmpi->data; 1416 c->ptap = ptap; 1417 ptap->api = NULL; 1418 ptap->apj = NULL; 1419 ptap->merge = merge; 1420 ptap->rmax = rmax; 1421 1422 *C = Cmpi; 1423 #if defined(PETSC_USE_INFO) 1424 if (bi[pn] != 0) { 1425 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1426 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1427 } else { 1428 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1429 } 1430 #endif 1431 PetscFunctionReturn(0); 1432 } 1433 1434 #undef __FUNCT__ 1435 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ" 1436 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1437 { 1438 PetscErrorCode ierr; 1439 Mat_Merge_SeqsToMPI *merge; 1440 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1441 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1442 Mat_PtAPMPI *ptap; 1443 PetscInt *adj; 1444 PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1445 MatScalar *ada,*ca,valtmp; 1446 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1447 MPI_Comm comm; 1448 PetscMPIInt size,rank,taga,*len_s; 1449 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1450 PetscInt **buf_ri,**buf_rj; 1451 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1452 MPI_Request *s_waits,*r_waits; 1453 MPI_Status *status; 1454 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1455 PetscInt *ai,*aj,*coi,*coj; 1456 PetscInt *poJ,*pdJ; 1457 Mat A_loc; 1458 Mat_SeqAIJ *a_loc; 1459 1460 PetscFunctionBegin; 1461 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1462 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1463 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1464 1465 ptap = c->ptap; 1466 merge = ptap->merge; 1467 1468 /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1469 /*------------------------------------------*/ 1470 /* get data from symbolic products */ 1471 coi = merge->coi; coj = merge->coj; 1472 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1473 bi = merge->bi; bj = merge->bj; 1474 owners = merge->rowmap->range; 1475 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1476 1477 /* get A_loc by taking all local rows of A */ 1478 A_loc = ptap->A_loc; 1479 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1480 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1481 ai = a_loc->i; 1482 aj = a_loc->j; 1483 1484 for (i=0; i<am; i++) { 1485 anz = ai[i+1] - ai[i]; 1486 adj = aj + ai[i]; 1487 ada = a_loc->a + ai[i]; 1488 1489 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1490 /*-------------------------------------------------------------*/ 1491 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1492 pnz = po->i[i+1] - po->i[i]; 1493 poJ = po->j + po->i[i]; 1494 pA = po->a + po->i[i]; 1495 for (j=0; j<pnz; j++) { 1496 row = poJ[j]; 1497 cj = coj + coi[row]; 1498 ca = coa + coi[row]; 1499 /* perform sparse axpy */ 1500 nexta = 0; 1501 valtmp = pA[j]; 1502 for (k=0; nexta<anz; k++) { 1503 if (cj[k] == adj[nexta]) { 1504 ca[k] += valtmp*ada[nexta]; 1505 nexta++; 1506 } 1507 } 1508 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1509 } 1510 1511 /* put the value into Cd (diagonal part) */ 1512 pnz = pd->i[i+1] - pd->i[i]; 1513 pdJ = pd->j + pd->i[i]; 1514 pA = pd->a + pd->i[i]; 1515 for (j=0; j<pnz; j++) { 1516 row = pdJ[j]; 1517 cj = bj + bi[row]; 1518 ca = ba + bi[row]; 1519 /* perform sparse axpy */ 1520 nexta = 0; 1521 valtmp = pA[j]; 1522 for (k=0; nexta<anz; k++) { 1523 if (cj[k] == adj[nexta]) { 1524 ca[k] += valtmp*ada[nexta]; 1525 nexta++; 1526 } 1527 } 1528 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1529 } 1530 } 1531 1532 /* 3) send and recv matrix values coa */ 1533 /*------------------------------------*/ 1534 buf_ri = merge->buf_ri; 1535 buf_rj = merge->buf_rj; 1536 len_s = merge->len_s; 1537 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1538 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1539 1540 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1541 for (proc=0,k=0; proc<size; proc++) { 1542 if (!len_s[proc]) continue; 1543 i = merge->owners_co[proc]; 1544 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1545 k++; 1546 } 1547 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1548 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1549 1550 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1551 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1552 ierr = PetscFree(coa);CHKERRQ(ierr); 1553 1554 /* 4) insert local Cseq and received values into Cmpi */ 1555 /*----------------------------------------------------*/ 1556 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1557 for (k=0; k<merge->nrecv; k++) { 1558 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1559 nrows = *(buf_ri_k[k]); 1560 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1561 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1562 } 1563 1564 for (i=0; i<cm; i++) { 1565 row = owners[rank] + i; /* global row index of C_seq */ 1566 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1567 ba_i = ba + bi[i]; 1568 bnz = bi[i+1] - bi[i]; 1569 /* add received vals into ba */ 1570 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1571 /* i-th row */ 1572 if (i == *nextrow[k]) { 1573 cnz = *(nextci[k]+1) - *nextci[k]; 1574 cj = buf_rj[k] + *(nextci[k]); 1575 ca = abuf_r[k] + *(nextci[k]); 1576 nextcj = 0; 1577 for (j=0; nextcj<cnz; j++) { 1578 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1579 ba_i[j] += ca[nextcj++]; 1580 } 1581 } 1582 nextrow[k]++; nextci[k]++; 1583 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1584 } 1585 } 1586 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1587 } 1588 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1589 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1590 1591 ierr = PetscFree(ba);CHKERRQ(ierr); 1592 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1593 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1594 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1595 PetscFunctionReturn(0); 1596 } 1597 1598 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*); 1599 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ(); 1600 differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */ 1601 #undef __FUNCT__ 1602 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ" 1603 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1604 { 1605 PetscErrorCode ierr; 1606 Mat Cmpi,A_loc,POt,PDt; 1607 Mat_PtAPMPI *ptap; 1608 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1609 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c; 1610 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1611 PetscInt nnz; 1612 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1613 PetscInt am =A->rmap->n,pn=P->cmap->n; 1614 MPI_Comm comm; 1615 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1616 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1617 PetscInt len,proc,*dnz,*onz,*owners; 1618 PetscInt nzi,*bi,*bj; 1619 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1620 MPI_Request *swaits,*rwaits; 1621 MPI_Status *sstatus,rstatus; 1622 Mat_Merge_SeqsToMPI *merge; 1623 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1624 PetscReal afill =1.0,afill_tmp; 1625 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax; 1626 PetscScalar *vals; 1627 Mat_SeqAIJ *a_loc,*pdt,*pot; 1628 PetscTable ta; 1629 1630 PetscFunctionBegin; 1631 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1632 /* check if matrix local sizes are compatible */ 1633 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 1634 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); 1635 } 1636 1637 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1638 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1639 1640 /* create struct Mat_PtAPMPI and attached it to C later */ 1641 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1642 1643 /* get A_loc by taking all local rows of A */ 1644 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1645 1646 ptap->A_loc = A_loc; 1647 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1648 ai = a_loc->i; 1649 aj = a_loc->j; 1650 1651 /* determine symbolic Co=(p->B)^T*A - send to others */ 1652 /*----------------------------------------------------*/ 1653 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1654 pdt = (Mat_SeqAIJ*)PDt->data; 1655 pdti = pdt->i; pdtj = pdt->j; 1656 1657 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1658 pot = (Mat_SeqAIJ*)POt->data; 1659 poti = pot->i; potj = pot->j; 1660 1661 /* then, compute symbolic Co = (p->B)^T*A */ 1662 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1663 >= (num of nonzero rows of C_seq) - pn */ 1664 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1665 coi[0] = 0; 1666 1667 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1668 nnz = fill*(poti[pon] + ai[am]); 1669 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1670 current_space = free_space; 1671 1672 /* create and initialize a linked list */ 1673 ierr = PetscTableCreate(2*a_loc->rmax,aN,&ta);CHKERRQ(ierr); 1674 MatRowMergeMax_SeqAIJ(a_loc,am,ta); 1675 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 1676 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 1677 1678 for (i=0; i<pon; i++) { 1679 pnz = poti[i+1] - poti[i]; 1680 ptJ = potj + poti[i]; 1681 for (j=0; j<pnz; j++) { 1682 row = ptJ[j]; /* row of A_loc == col of Pot */ 1683 anz = ai[row+1] - ai[row]; 1684 Jptr = aj + ai[row]; 1685 /* add non-zero cols of AP into the sorted linked list lnk */ 1686 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1687 } 1688 nnz = lnk[0]; 1689 1690 /* If free space is not available, double the total space in the list */ 1691 if (current_space->local_remaining<nnz) { 1692 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1693 nspacedouble++; 1694 } 1695 1696 /* Copy data into free space, and zero out denserows */ 1697 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1698 1699 current_space->array += nnz; 1700 current_space->local_used += nnz; 1701 current_space->local_remaining -= nnz; 1702 1703 coi[i+1] = coi[i] + nnz; 1704 } 1705 1706 ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1707 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1708 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */ 1709 1710 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1711 if (afill_tmp > afill) afill = afill_tmp; 1712 1713 /* send j-array (coj) of Co to other processors */ 1714 /*----------------------------------------------*/ 1715 /* determine row ownership */ 1716 ierr = PetscNew(&merge);CHKERRQ(ierr); 1717 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1718 1719 merge->rowmap->n = pn; 1720 merge->rowmap->bs = 1; 1721 1722 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1723 owners = merge->rowmap->range; 1724 1725 /* determine the number of messages to send, their lengths */ 1726 ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1727 ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 1728 1729 len_s = merge->len_s; 1730 merge->nsend = 0; 1731 1732 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1733 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1734 1735 proc = 0; 1736 for (i=0; i<pon; i++) { 1737 while (prmap[i] >= owners[proc+1]) proc++; 1738 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1739 len_s[proc] += coi[i+1] - coi[i]; 1740 } 1741 1742 len = 0; /* max length of buf_si[] */ 1743 owners_co[0] = 0; 1744 for (proc=0; proc<size; proc++) { 1745 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1746 if (len_si[proc]) { 1747 merge->nsend++; 1748 len_si[proc] = 2*(len_si[proc] + 1); 1749 len += len_si[proc]; 1750 } 1751 } 1752 1753 /* determine the number and length of messages to receive for coi and coj */ 1754 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1755 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1756 1757 /* post the Irecv and Isend of coj */ 1758 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1759 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1760 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1761 for (proc=0, k=0; proc<size; proc++) { 1762 if (!len_s[proc]) continue; 1763 i = owners_co[proc]; 1764 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1765 k++; 1766 } 1767 1768 /* receives and sends of coj are complete */ 1769 ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1770 for (i=0; i<merge->nrecv; i++) { 1771 PetscMPIInt icompleted; 1772 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1773 } 1774 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1775 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1776 1777 /* add received column indices into table to update Armax */ 1778 for (k=0; k<merge->nrecv; k++) {/* k-th received message */ 1779 Jptr = buf_rj[k]; 1780 for (j=0; j<merge->len_r[k]; j++) { 1781 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1782 } 1783 } 1784 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 1785 1786 /* send and recv coi */ 1787 /*-------------------*/ 1788 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1789 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1790 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1791 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1792 for (proc=0,k=0; proc<size; proc++) { 1793 if (!len_s[proc]) continue; 1794 /* form outgoing message for i-structure: 1795 buf_si[0]: nrows to be sent 1796 [1:nrows]: row index (global) 1797 [nrows+1:2*nrows+1]: i-structure index 1798 */ 1799 /*-------------------------------------------*/ 1800 nrows = len_si[proc]/2 - 1; 1801 buf_si_i = buf_si + nrows+1; 1802 buf_si[0] = nrows; 1803 buf_si_i[0] = 0; 1804 nrows = 0; 1805 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1806 nzi = coi[i+1] - coi[i]; 1807 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1808 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1809 nrows++; 1810 } 1811 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1812 k++; 1813 buf_si += len_si[proc]; 1814 } 1815 i = merge->nrecv; 1816 while (i--) { 1817 PetscMPIInt icompleted; 1818 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1819 } 1820 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1821 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1822 ierr = PetscFree(len_si);CHKERRQ(ierr); 1823 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1824 ierr = PetscFree(swaits);CHKERRQ(ierr); 1825 ierr = PetscFree(sstatus);CHKERRQ(ierr); 1826 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1827 1828 /* compute the local portion of C (mpi mat) */ 1829 /*------------------------------------------*/ 1830 /* allocate bi array and free space for accumulating nonzero column info */ 1831 ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1832 bi[0] = 0; 1833 1834 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1835 nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1836 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1837 current_space = free_space; 1838 1839 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1840 for (k=0; k<merge->nrecv; k++) { 1841 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1842 nrows = *buf_ri_k[k]; 1843 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1844 nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure */ 1845 } 1846 1847 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 1848 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1849 rmax = 0; 1850 for (i=0; i<pn; i++) { 1851 /* add pdt[i,:]*AP into lnk */ 1852 pnz = pdti[i+1] - pdti[i]; 1853 ptJ = pdtj + pdti[i]; 1854 for (j=0; j<pnz; j++) { 1855 row = ptJ[j]; /* row of AP == col of Pt */ 1856 anz = ai[row+1] - ai[row]; 1857 Jptr = aj + ai[row]; 1858 /* add non-zero cols of AP into the sorted linked list lnk */ 1859 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1860 } 1861 1862 /* add received col data into lnk */ 1863 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1864 if (i == *nextrow[k]) { /* i-th row */ 1865 nzi = *(nextci[k]+1) - *nextci[k]; 1866 Jptr = buf_rj[k] + *nextci[k]; 1867 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 1868 nextrow[k]++; nextci[k]++; 1869 } 1870 } 1871 nnz = lnk[0]; 1872 1873 /* if free space is not available, make more free space */ 1874 if (current_space->local_remaining<nnz) { 1875 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1876 nspacedouble++; 1877 } 1878 /* copy data into free space, then initialize lnk */ 1879 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1880 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1881 1882 current_space->array += nnz; 1883 current_space->local_used += nnz; 1884 current_space->local_remaining -= nnz; 1885 1886 bi[i+1] = bi[i] + nnz; 1887 if (nnz > rmax) rmax = nnz; 1888 } 1889 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1890 1891 ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 1892 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1893 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 1894 if (afill_tmp > afill) afill = afill_tmp; 1895 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 1896 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 1897 1898 ierr = MatDestroy(&POt);CHKERRQ(ierr); 1899 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1900 1901 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 1902 /*----------------------------------------------------------------------------------*/ 1903 ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 1904 1905 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1906 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1907 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 1908 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1909 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1910 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1911 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1912 for (i=0; i<pn; i++) { 1913 row = i + rstart; 1914 nnz = bi[i+1] - bi[i]; 1915 Jptr = bj + bi[i]; 1916 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1917 } 1918 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1919 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1920 ierr = PetscFree(vals);CHKERRQ(ierr); 1921 1922 merge->bi = bi; 1923 merge->bj = bj; 1924 merge->coi = coi; 1925 merge->coj = coj; 1926 merge->buf_ri = buf_ri; 1927 merge->buf_rj = buf_rj; 1928 merge->owners_co = owners_co; 1929 merge->destroy = Cmpi->ops->destroy; 1930 merge->duplicate = Cmpi->ops->duplicate; 1931 1932 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 1933 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1934 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1935 1936 /* attach the supporting struct to Cmpi for reuse */ 1937 c = (Mat_MPIAIJ*)Cmpi->data; 1938 1939 c->ptap = ptap; 1940 ptap->api = NULL; 1941 ptap->apj = NULL; 1942 ptap->merge = merge; 1943 ptap->rmax = rmax; 1944 ptap->apa = NULL; 1945 1946 *C = Cmpi; 1947 #if defined(PETSC_USE_INFO) 1948 if (bi[pn] != 0) { 1949 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1950 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1951 } else { 1952 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1953 } 1954 #endif 1955 PetscFunctionReturn(0); 1956 } 1957