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