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 #include <petsc/private/vecscatterimpl.h> 13 14 #if defined(PETSC_HAVE_HYPRE) 15 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 16 #endif 17 18 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 19 { 20 PetscErrorCode ierr; 21 #if defined(PETSC_HAVE_HYPRE) 22 const char *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"}; 23 PetscInt nalg = 4; 24 #else 25 const char *algTypes[3] = {"scalable","nonscalable","seqmpi"}; 26 PetscInt nalg = 3; 27 #endif 28 PetscInt alg = 1; /* set nonscalable algorithm as default */ 29 MPI_Comm comm; 30 PetscBool flg; 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 = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatMatMult","Mat");CHKERRQ(ierr); 38 ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr); 39 ierr = PetscOptionsEnd();CHKERRQ(ierr); 40 41 if (!flg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */ 42 MatInfo Ainfo,Binfo; 43 PetscInt nz_local; 44 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 45 46 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 47 ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr); 48 nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated); 49 50 if (B->cmap->N > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 51 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 52 53 if (alg_scalable) { 54 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 55 ierr = PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,fill*nz_local);CHKERRQ(ierr); 56 } 57 } 58 59 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 60 switch (alg) { 61 case 1: 62 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr); 63 break; 64 case 2: 65 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);CHKERRQ(ierr); 66 break; 67 #if defined(PETSC_HAVE_HYPRE) 68 case 3: 69 ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr); 70 break; 71 #endif 72 default: 73 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 74 break; 75 } 76 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 77 78 if (alg == 0 || alg == 1) { 79 Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 80 Mat_APMPI *ap = c->ap; 81 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 82 ap->freestruct = PETSC_FALSE; 83 ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 84 ierr = PetscOptionsEnd();CHKERRQ(ierr); 85 } 86 } 87 88 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 89 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 90 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 95 { 96 PetscErrorCode ierr; 97 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 98 Mat_APMPI *ptap = a->ap; 99 100 PetscFunctionBegin; 101 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 102 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 103 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 104 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 105 ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr); 106 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 107 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 108 ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 109 ierr = ptap->destroy(A);CHKERRQ(ierr); 110 ierr = PetscFree(ptap);CHKERRQ(ierr); 111 PetscFunctionReturn(0); 112 } 113 114 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C) 115 { 116 PetscErrorCode ierr; 117 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 118 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 119 Mat_SeqAIJ *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 120 PetscScalar *cda=cd->a,*coa=co->a; 121 Mat_SeqAIJ *p_loc,*p_oth; 122 PetscScalar *apa,*ca; 123 PetscInt cm =C->rmap->n; 124 Mat_APMPI *ptap=c->ap; 125 PetscInt *api,*apj,*apJ,i,k; 126 PetscInt cstart=C->cmap->rstart; 127 PetscInt cdnz,conz,k0,k1; 128 MPI_Comm comm; 129 PetscMPIInt size; 130 131 PetscFunctionBegin; 132 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 133 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 134 135 if (!ptap->P_oth && size>1) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 136 137 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 138 /*-----------------------------------------------------*/ 139 /* update numerical values of P_oth and P_loc */ 140 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 141 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 142 143 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 144 /*----------------------------------------------------------*/ 145 /* get data from symbolic products */ 146 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 147 p_oth = NULL; 148 if (size >1) { 149 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 150 } 151 152 /* get apa for storing dense row A[i,:]*P */ 153 apa = ptap->apa; 154 155 api = ptap->api; 156 apj = ptap->apj; 157 for (i=0; i<cm; i++) { 158 /* compute apa = A[i,:]*P */ 159 AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 160 161 /* set values in C */ 162 apJ = apj + api[i]; 163 cdnz = cd->i[i+1] - cd->i[i]; 164 conz = co->i[i+1] - co->i[i]; 165 166 /* 1st off-diagoanl part of C */ 167 ca = coa + co->i[i]; 168 k = 0; 169 for (k0=0; k0<conz; k0++) { 170 if (apJ[k] >= cstart) break; 171 ca[k0] = apa[apJ[k]]; 172 apa[apJ[k++]] = 0.0; 173 } 174 175 /* diagonal part of C */ 176 ca = cda + cd->i[i]; 177 for (k1=0; k1<cdnz; k1++) { 178 ca[k1] = apa[apJ[k]]; 179 apa[apJ[k++]] = 0.0; 180 } 181 182 /* 2nd off-diagoanl part of C */ 183 ca = coa + co->i[i]; 184 for (; k0<conz; k0++) { 185 ca[k0] = apa[apJ[k]]; 186 apa[apJ[k++]] = 0.0; 187 } 188 } 189 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 190 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 191 192 if (ptap->freestruct) { 193 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 194 } 195 PetscFunctionReturn(0); 196 } 197 198 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C) 199 { 200 PetscErrorCode ierr; 201 MPI_Comm comm; 202 PetscMPIInt size; 203 Mat Cmpi; 204 Mat_APMPI *ptap; 205 PetscFreeSpaceList free_space=NULL,current_space=NULL; 206 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c; 207 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 208 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 209 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 210 PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 211 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 212 PetscBT lnkbt; 213 PetscReal afill; 214 MatType mtype; 215 216 PetscFunctionBegin; 217 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 218 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 219 220 /* create struct Mat_APMPI and attached it to C later */ 221 ierr = PetscNew(&ptap);CHKERRQ(ierr); 222 223 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 224 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 225 226 /* get P_loc by taking all local rows of P */ 227 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 228 229 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 230 pi_loc = p_loc->i; pj_loc = p_loc->j; 231 if (size > 1) { 232 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 233 pi_oth = p_oth->i; pj_oth = p_oth->j; 234 } else { 235 p_oth = NULL; 236 pi_oth = NULL; pj_oth = NULL; 237 } 238 239 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 240 /*-------------------------------------------------------------------*/ 241 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 242 ptap->api = api; 243 api[0] = 0; 244 245 /* create and initialize a linked list */ 246 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 247 248 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 249 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 250 current_space = free_space; 251 252 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 253 for (i=0; i<am; i++) { 254 /* diagonal portion of A */ 255 nzi = adi[i+1] - adi[i]; 256 for (j=0; j<nzi; j++) { 257 row = *adj++; 258 pnz = pi_loc[row+1] - pi_loc[row]; 259 Jptr = pj_loc + pi_loc[row]; 260 /* add non-zero cols of P into the sorted linked list lnk */ 261 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 262 } 263 /* off-diagonal portion of A */ 264 nzi = aoi[i+1] - aoi[i]; 265 for (j=0; j<nzi; j++) { 266 row = *aoj++; 267 pnz = pi_oth[row+1] - pi_oth[row]; 268 Jptr = pj_oth + pi_oth[row]; 269 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 270 } 271 272 apnz = lnk[0]; 273 api[i+1] = api[i] + apnz; 274 275 /* if free space is not available, double the total space in the list */ 276 if (current_space->local_remaining<apnz) { 277 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 278 nspacedouble++; 279 } 280 281 /* Copy data into free space, then initialize lnk */ 282 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 283 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 284 285 current_space->array += apnz; 286 current_space->local_used += apnz; 287 current_space->local_remaining -= apnz; 288 } 289 290 /* Allocate space for apj, initialize apj, and */ 291 /* destroy list of free space and other temporary array(s) */ 292 ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 293 apj = ptap->apj; 294 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 295 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 296 297 /* malloc apa to store dense row A[i,:]*P */ 298 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 299 300 /* create and assemble symbolic parallel matrix Cmpi */ 301 /*----------------------------------------------------*/ 302 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 303 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 304 ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 305 306 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 307 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 308 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 309 310 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 311 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 312 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 313 ierr = MatPreallocateFinalize(dnz,onz);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->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 320 321 /* attach the supporting struct to Cmpi for reuse */ 322 c = (Mat_MPIAIJ*)Cmpi->data; 323 c->ap = 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 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 346 { 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 if (scall == MAT_INITIAL_MATRIX) { 351 *C = NULL; 352 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 353 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 354 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 355 } 356 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 357 ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 358 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 359 PetscFunctionReturn(0); 360 } 361 362 typedef struct { 363 Mat workB,Bb,Cb,workB1,Bb1,Cb1; 364 MPI_Request *rwaits,*swaits; 365 PetscInt numBb; /* num of Bb matrices */ 366 PetscInt nsends,nrecvs; 367 MPI_Datatype *stype,*rtype; 368 } MPIAIJ_MPIDense; 369 370 PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx) 371 { 372 MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*)ctx; 373 PetscErrorCode ierr; 374 PetscInt i; 375 376 PetscFunctionBegin; 377 ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 378 379 if (contents->numBb) { 380 ierr = MatDestroy(&contents->Bb);CHKERRQ(ierr); 381 ierr = MatDestroy(&contents->Cb);CHKERRQ(ierr); 382 383 ierr = MatDestroy(&contents->workB1);CHKERRQ(ierr); 384 ierr = MatDestroy(&contents->Bb1);CHKERRQ(ierr); 385 ierr = MatDestroy(&contents->Cb1);CHKERRQ(ierr); 386 } 387 for (i=0; i<contents->nsends; i++) { 388 ierr = MPI_Type_free(&contents->stype[i]);CHKERRQ(ierr); 389 } 390 for (i=0; i<contents->nrecvs; i++) { 391 ierr = MPI_Type_free(&contents->rtype[i]);CHKERRQ(ierr); 392 } 393 ierr = PetscFree4(contents->stype,contents->rtype,contents->rwaits,contents->swaits);CHKERRQ(ierr); 394 ierr = PetscFree(contents);CHKERRQ(ierr); 395 PetscFunctionReturn(0); 396 } 397 398 /* 399 This is a "dummy function" that handles the case where matrix C was created as a dense matrix 400 directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option 401 402 It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C 403 */ 404 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C) 405 { 406 PetscBool flg; 407 PetscErrorCode ierr; 408 409 PetscFunctionBegin; 410 ierr = PetscObjectTypeCompare((PetscObject)A,MATNEST,&flg);CHKERRQ(ierr); 411 if (flg) { 412 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 413 ierr = MatMatMultSymbolic_Nest_Dense(A,B,PETSC_DEFAULT,&C);CHKERRQ(ierr); 414 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 415 C->ops->matmultnumeric = MatMatMultNumeric_Nest_Dense; 416 } else { 417 ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,PETSC_DEFAULT,&C);CHKERRQ(ierr); 418 } 419 ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr); 420 PetscFunctionReturn(0); 421 } 422 423 /* 424 Create Bb, Cb, Bb1 and Cb1 matrices to be used by MatMatMult_MPIAIJ_MPIDense(). 425 These matrices are used as wrappers for sub-columns of B and C, thus their own matrix operations are not used. 426 Modified from MatCreateDense(). 427 */ 428 PETSC_STATIC_INLINE PetscErrorCode MatCreateSubMPIDense_private(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt rbs,PetscInt cbs,PetscScalar *data,Mat *A) 429 { 430 PetscErrorCode ierr; 431 432 PetscFunctionBegin; 433 ierr = MatCreate(comm,A);CHKERRQ(ierr); 434 ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 435 ierr = MatSetBlockSizes(*A,rbs,cbs);CHKERRQ(ierr); 436 ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 437 ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 438 (*A)->assembled = PETSC_TRUE; 439 PetscFunctionReturn(0); 440 } 441 442 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 443 { 444 PetscErrorCode ierr; 445 Mat_MPIAIJ *aij=(Mat_MPIAIJ*)A->data; 446 Mat_MPIDense *b=(Mat_MPIDense*)B->data; 447 Mat_SeqDense *bseq=(Mat_SeqDense*)(b->A)->data; 448 PetscInt nz=aij->B->cmap->n,nsends,nrecvs,i,nrows_to,j,lda=bseq->lda; 449 PetscContainer container; 450 MPIAIJ_MPIDense *contents; 451 VecScatter ctx=aij->Mvctx; 452 PetscInt Am=A->rmap->n,Bm=B->rmap->n,Bn=B->cmap->n,BN=B->cmap->N,Bbn,Bbn1,bs,nrows_from; 453 MPI_Comm comm; 454 MPI_Datatype type1,*stype,*rtype; 455 const PetscInt *sindices,*sstarts,*rstarts; 456 PetscMPIInt *disp; 457 458 PetscFunctionBegin; 459 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 460 if (!(*C)) { 461 ierr = MatCreate(comm,C);CHKERRQ(ierr); 462 ierr = MatSetSizes(*C,Am,Bn,A->rmap->N,BN);CHKERRQ(ierr); 463 ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 464 ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 465 ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr); 466 ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 467 ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 468 } else { 469 /* Check matrix size */ 470 if ((*C)->rmap->n != A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local row dimensions are incompatible, %D != %D",(*C)->rmap->n,A->rmap->n); 471 if ((*C)->cmap->n != B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local column dimensions are incompatible, %D != %D",(*C)->cmap->n,B->cmap->n); 472 } 473 474 (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 475 476 ierr = PetscNew(&contents);CHKERRQ(ierr); 477 contents->numBb = 0; 478 479 ierr = VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);CHKERRQ(ierr); 480 ierr = VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);CHKERRQ(ierr); 481 482 /* Create column block of B and C for memory scalability when BN is too large */ 483 /* Estimate Bbn, column size of Bb */ 484 if (nz) { 485 Bbn1 = 2*Am*BN/nz; 486 } else Bbn1 = BN; 487 488 bs = PetscAbs(B->cmap->bs); 489 Bbn1 = Bbn1/bs *bs; /* Bbn1 is a multiple of bs */ 490 if (Bbn1 > BN) Bbn1 = BN; 491 ierr = MPI_Allreduce(&Bbn1,&Bbn,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 492 493 /* Enable runtime option for Bbn */ 494 ierr = PetscOptionsBegin(comm,((PetscObject)*C)->prefix,"MatMatMult","Mat");CHKERRQ(ierr); 495 ierr = PetscOptionsInt("-matmatmult_Bbn","Number of columns in Bb","MatMatMult",Bbn,&Bbn,NULL);CHKERRQ(ierr); 496 if (Bbn > BN) SETERRQ2(comm,PETSC_ERR_ARG_SIZ,"Bbn=%D cannot be larger than %D, column size of B",Bbn,BN); 497 ierr = PetscOptionsEnd();CHKERRQ(ierr); 498 499 if (Bbn < BN) { 500 contents->numBb = BN/Bbn; 501 Bbn1 = BN - contents->numBb*Bbn; 502 } 503 504 if (contents->numBb) { 505 PetscScalar data[1]; /* fake array for Bb and Cb */ 506 ierr = PetscInfo3(*C,"use Bb, BN=%D, Bbn=%D; numBb=%D\n",BN,Bbn,contents->numBb);CHKERRQ(ierr); 507 ierr = MatCreateSubMPIDense_private(comm,B->rmap->n,PETSC_DECIDE,A->rmap->N,Bbn,B->rmap->bs,B->cmap->bs,data,&contents->Bb);CHKERRQ(ierr); 508 ierr = MatCreateSubMPIDense_private(comm,Am,PETSC_DECIDE,A->rmap->N,Bbn,(*C)->rmap->bs,(*C)->cmap->bs,data,&contents->Cb);CHKERRQ(ierr); 509 510 if (Bbn1) { /* Create Bb1 and Cb1 for the remaining columns */ 511 ierr = PetscInfo2(*C,"use Bb1, BN=%D, Bbn1=%D\n",BN,Bbn1);CHKERRQ(ierr); 512 ierr = MatCreateSubMPIDense_private(comm,B->rmap->n,PETSC_DECIDE,A->rmap->N,Bbn1,B->rmap->bs,B->cmap->bs,data,&contents->Bb1);CHKERRQ(ierr); 513 ierr = MatCreateSubMPIDense_private(comm,Am,PETSC_DECIDE,A->rmap->N,Bbn1,(*C)->rmap->bs,(*C)->cmap->bs,data,&contents->Cb1);CHKERRQ(ierr); 514 515 /* Create work matrix used to store off processor rows of B needed for local product */ 516 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn1,NULL,&contents->workB1);CHKERRQ(ierr); 517 } 518 } 519 520 /* Create work matrix used to store off processor rows of B needed for local product */ 521 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn,NULL,&contents->workB);CHKERRQ(ierr); 522 523 /* Use MPI derived data type to reduce memory required by the send/recv buffers */ 524 ierr = PetscMalloc4(nsends,&stype,nrecvs,&rtype,nrecvs,&contents->rwaits,nsends,&contents->swaits);CHKERRQ(ierr); 525 contents->stype = stype; 526 contents->nsends = nsends; 527 528 contents->rtype = rtype; 529 contents->nrecvs = nrecvs; 530 531 ierr = PetscMalloc1(Bm+1,&disp);CHKERRQ(ierr); 532 for (i=0; i<nsends; i++) { 533 nrows_to = sstarts[i+1]-sstarts[i]; 534 for (j=0; j<nrows_to; j++){ 535 disp[j] = sindices[sstarts[i]+j]; /* rowB to be sent */ 536 } 537 ierr = MPI_Type_create_indexed_block(nrows_to,1,(const PetscMPIInt *)disp,MPIU_SCALAR,&type1);CHKERRQ(ierr); 538 539 ierr = MPI_Type_create_resized(type1,0,lda*sizeof(PetscScalar),&stype[i]);CHKERRQ(ierr); 540 ierr = MPI_Type_commit(&stype[i]);CHKERRQ(ierr); 541 ierr = MPI_Type_free(&type1);CHKERRQ(ierr); 542 } 543 544 for (i=0; i<nrecvs; i++) { 545 /* received values from a process form a (nrows_from x Bbn) row block in workB (column-wise) */ 546 nrows_from = rstarts[i+1]-rstarts[i]; 547 disp[0] = 0; 548 ierr = MPI_Type_create_indexed_block(1, nrows_from, (const PetscMPIInt *)disp, MPIU_SCALAR, &type1);CHKERRQ(ierr); 549 ierr = MPI_Type_create_resized(type1, 0, nz*sizeof(PetscScalar), &rtype[i]);CHKERRQ(ierr); 550 ierr = MPI_Type_commit(&rtype[i]);CHKERRQ(ierr); 551 ierr = MPI_Type_free(&type1);CHKERRQ(ierr); 552 } 553 554 ierr = PetscFree(disp);CHKERRQ(ierr); 555 ierr = VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);CHKERRQ(ierr); 556 ierr = VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);CHKERRQ(ierr); 557 558 ierr = PetscContainerCreate(comm,&container);CHKERRQ(ierr); 559 ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 560 ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 561 ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 562 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 563 PetscFunctionReturn(0); 564 } 565 566 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 567 /* 568 Performs an efficient scatter on the rows of B needed by this process; this is 569 a modification of the VecScatterBegin_() routines. 570 571 Input: Bbidx = 0: B = Bb 572 = 1: B = Bb1, see MatMatMultSymbolic_MPIAIJ_MPIDense() 573 */ 574 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,PetscInt Bbidx,Mat C,Mat *outworkB) 575 { 576 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 577 PetscErrorCode ierr; 578 const PetscScalar *b; 579 PetscScalar *rvalues; 580 VecScatter ctx = aij->Mvctx; 581 const PetscInt *sindices,*sstarts,*rstarts; 582 const PetscMPIInt *sprocs,*rprocs; 583 PetscInt i,nsends,nrecvs,nrecvs2; 584 MPI_Request *swaits,*rwaits; 585 MPI_Comm comm; 586 PetscMPIInt tag=((PetscObject)ctx)->tag,ncols=B->cmap->N,nrows=aij->B->cmap->n,imdex,nsends_mpi,nrecvs_mpi; 587 MPI_Status status; 588 MPIAIJ_MPIDense *contents; 589 PetscContainer container; 590 Mat workB; 591 MPI_Datatype *stype,*rtype; 592 593 PetscFunctionBegin; 594 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 595 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 596 if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 597 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 598 599 ierr = VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);CHKERRQ(ierr); 600 ierr = VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL/*bs*/);CHKERRQ(ierr); 601 ierr = PetscMPIIntCast(nsends,&nsends_mpi);CHKERRQ(ierr); 602 ierr = PetscMPIIntCast(nrecvs,&nrecvs_mpi);CHKERRQ(ierr); 603 if (Bbidx == 0) { 604 workB = *outworkB = contents->workB; 605 } else { 606 workB = *outworkB = contents->workB1; 607 } 608 if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",workB->cmap->n,nrows); 609 swaits = contents->swaits; 610 rwaits = contents->rwaits; 611 612 ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr); 613 ierr = MatDenseGetArray(workB,&rvalues);CHKERRQ(ierr); 614 615 /* Post recv, use MPI derived data type to save memory */ 616 rtype = contents->rtype; 617 for (i=0; i<nrecvs; i++) { 618 ierr = MPI_Irecv(rvalues+(rstarts[i]-rstarts[0]),ncols,rtype[i],rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 619 } 620 621 stype = contents->stype; 622 for (i=0; i<nsends; i++) { 623 ierr = MPI_Isend(b,ncols,stype[i],sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 624 } 625 626 nrecvs2 = nrecvs; 627 while (nrecvs2) { 628 ierr = MPI_Waitany(nrecvs_mpi,rwaits,&imdex,&status);CHKERRQ(ierr); 629 nrecvs2--; 630 } 631 if (nsends) {ierr = MPI_Waitall(nsends_mpi,swaits,MPI_STATUSES_IGNORE);CHKERRQ(ierr);} 632 633 ierr = VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL);CHKERRQ(ierr); 634 ierr = VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL);CHKERRQ(ierr); 635 ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr); 636 ierr = MatDenseRestoreArray(workB,&rvalues);CHKERRQ(ierr); 637 ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 638 ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 639 PetscFunctionReturn(0); 640 } 641 642 /* 643 Compute Cb = A*Bb 644 */ 645 PETSC_STATIC_INLINE PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense_private(Mat A,Mat Bb,PetscInt Bbidx,PetscInt start,Mat C,const PetscScalar *barray,PetscScalar *carray,Mat Cb) 646 { 647 PetscErrorCode ierr; 648 PetscInt start1; 649 Mat workB; 650 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 651 Mat_MPIDense *cbdense = (Mat_MPIDense*)Cb->data; 652 653 PetscFunctionBegin; 654 /* Place barray to Bb */ 655 start1 = start*Bb->rmap->n; 656 ierr = MatDensePlaceArray(Bb,barray+start1);CHKERRQ(ierr); 657 658 /* get off processor parts of Bb needed to complete Cb=A*Bb */ 659 ierr = MatMPIDenseScatter(A,Bb,Bbidx,C,&workB);CHKERRQ(ierr); 660 ierr = MatDenseResetArray(Bb);CHKERRQ(ierr); 661 662 /* off-diagonal block of A times nonlocal rows of Bb */ 663 /* Place carray to Cb */ 664 start1 = start*Cb->rmap->n; 665 ierr = MatDensePlaceArray(Cb,carray+start1);CHKERRQ(ierr); 666 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cbdense->A);CHKERRQ(ierr); 667 ierr = MatDenseResetArray(Cb);CHKERRQ(ierr); 668 PetscFunctionReturn(0); 669 } 670 671 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 672 { 673 PetscErrorCode ierr; 674 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 675 Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 676 Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 677 Mat workB; 678 MPIAIJ_MPIDense *contents; 679 PetscContainer container; 680 MPI_Comm comm; 681 PetscInt numBb; 682 683 PetscFunctionBegin; 684 /* diagonal block of A times all local rows of B*/ 685 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 686 687 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 688 ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 689 if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 690 ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 691 numBb = contents->numBb; 692 693 if (!numBb) { 694 /* get off processor parts of B needed to complete C=A*B */ 695 ierr = MatMPIDenseScatter(A,B,0,C,&workB);CHKERRQ(ierr); 696 697 /* off-diagonal block of A times nonlocal rows of B */ 698 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 699 700 } else { 701 const PetscScalar *barray; 702 PetscScalar *carray; 703 Mat Bb=contents->Bb,Cb=contents->Cb; 704 PetscInt BbN=Bb->cmap->N,start,i; 705 706 ierr = MatDenseGetArrayRead(B,&barray);CHKERRQ(ierr); 707 ierr = MatDenseGetArray(C,&carray);CHKERRQ(ierr); 708 for (i=0; i<numBb; i++) { 709 start = i*BbN; 710 ierr = MatMatMultNumeric_MPIAIJ_MPIDense_private(A,Bb,0,start,C,barray,carray,Cb);CHKERRQ(ierr); 711 } 712 713 if (contents->Bb1) { 714 Bb = contents->Bb1; Cb = contents->Cb1; 715 start = i*BbN; 716 ierr = MatMatMultNumeric_MPIAIJ_MPIDense_private(A,Bb,1,start,C,barray,carray,Cb);CHKERRQ(ierr); 717 } 718 ierr = MatDenseRestoreArrayRead(B,&barray);CHKERRQ(ierr); 719 ierr = MatDenseRestoreArray(C,&carray);CHKERRQ(ierr); 720 } 721 PetscFunctionReturn(0); 722 } 723 724 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 725 { 726 PetscErrorCode ierr; 727 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 728 Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 729 Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 730 PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj; 731 PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 732 Mat_SeqAIJ *p_loc,*p_oth; 733 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 734 PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 735 PetscInt cm = C->rmap->n,anz,pnz; 736 Mat_APMPI *ptap = c->ap; 737 PetscScalar *apa_sparse; 738 PetscInt *api,*apj,*apJ,i,j,k,row; 739 PetscInt cstart = C->cmap->rstart; 740 PetscInt cdnz,conz,k0,k1,nextp; 741 MPI_Comm comm; 742 PetscMPIInt size; 743 744 PetscFunctionBegin; 745 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 746 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 747 748 if (!ptap->P_oth && size>1) { 749 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 750 } 751 apa_sparse = ptap->apa; 752 753 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 754 /*-----------------------------------------------------*/ 755 /* update numerical values of P_oth and P_loc */ 756 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 757 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 758 759 /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 760 /*----------------------------------------------------------*/ 761 /* get data from symbolic products */ 762 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 763 pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a; 764 if (size >1) { 765 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 766 pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a; 767 } else { 768 p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL; 769 } 770 771 api = ptap->api; 772 apj = ptap->apj; 773 for (i=0; i<cm; i++) { 774 apJ = apj + api[i]; 775 776 /* diagonal portion of A */ 777 anz = adi[i+1] - adi[i]; 778 adj = ad->j + adi[i]; 779 ada = ad->a + adi[i]; 780 for (j=0; j<anz; j++) { 781 row = adj[j]; 782 pnz = pi_loc[row+1] - pi_loc[row]; 783 pj = pj_loc + pi_loc[row]; 784 pa = pa_loc + pi_loc[row]; 785 /* perform sparse axpy */ 786 valtmp = ada[j]; 787 nextp = 0; 788 for (k=0; nextp<pnz; k++) { 789 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 790 apa_sparse[k] += valtmp*pa[nextp++]; 791 } 792 } 793 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 794 } 795 796 /* off-diagonal portion of A */ 797 anz = aoi[i+1] - aoi[i]; 798 aoj = ao->j + aoi[i]; 799 aoa = ao->a + aoi[i]; 800 for (j=0; j<anz; j++) { 801 row = aoj[j]; 802 pnz = pi_oth[row+1] - pi_oth[row]; 803 pj = pj_oth + pi_oth[row]; 804 pa = pa_oth + pi_oth[row]; 805 /* perform sparse axpy */ 806 valtmp = aoa[j]; 807 nextp = 0; 808 for (k=0; nextp<pnz; k++) { 809 if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 810 apa_sparse[k] += valtmp*pa[nextp++]; 811 } 812 } 813 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 814 } 815 816 /* set values in C */ 817 cdnz = cd->i[i+1] - cd->i[i]; 818 conz = co->i[i+1] - co->i[i]; 819 820 /* 1st off-diagoanl part of C */ 821 ca = coa + co->i[i]; 822 k = 0; 823 for (k0=0; k0<conz; k0++) { 824 if (apJ[k] >= cstart) break; 825 ca[k0] = apa_sparse[k]; 826 apa_sparse[k] = 0.0; 827 k++; 828 } 829 830 /* diagonal part of C */ 831 ca = cda + cd->i[i]; 832 for (k1=0; k1<cdnz; k1++) { 833 ca[k1] = apa_sparse[k]; 834 apa_sparse[k] = 0.0; 835 k++; 836 } 837 838 /* 2nd off-diagoanl part of C */ 839 ca = coa + co->i[i]; 840 for (; k0<conz; k0++) { 841 ca[k0] = apa_sparse[k]; 842 apa_sparse[k] = 0.0; 843 k++; 844 } 845 } 846 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 847 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 848 849 if (ptap->freestruct) { 850 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 851 } 852 PetscFunctionReturn(0); 853 } 854 855 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */ 856 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 857 { 858 PetscErrorCode ierr; 859 MPI_Comm comm; 860 PetscMPIInt size; 861 Mat Cmpi; 862 Mat_APMPI *ptap; 863 PetscFreeSpaceList free_space = NULL,current_space=NULL; 864 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c; 865 Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 866 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 867 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 868 PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 869 PetscInt am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20; 870 PetscReal afill; 871 MatType mtype; 872 873 PetscFunctionBegin; 874 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 875 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 876 877 /* create struct Mat_APMPI and attached it to C later */ 878 ierr = PetscNew(&ptap);CHKERRQ(ierr); 879 880 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 881 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 882 883 /* get P_loc by taking all local rows of P */ 884 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 885 886 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 887 pi_loc = p_loc->i; pj_loc = p_loc->j; 888 if (size > 1) { 889 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 890 pi_oth = p_oth->i; pj_oth = p_oth->j; 891 } else { 892 p_oth = NULL; 893 pi_oth = NULL; pj_oth = NULL; 894 } 895 896 /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 897 /*-------------------------------------------------------------------*/ 898 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 899 ptap->api = api; 900 api[0] = 0; 901 902 ierr = PetscLLCondensedCreate_Scalable(lsize,&lnk);CHKERRQ(ierr); 903 904 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 905 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 906 current_space = free_space; 907 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 908 for (i=0; i<am; i++) { 909 /* diagonal portion of A */ 910 nzi = adi[i+1] - adi[i]; 911 for (j=0; j<nzi; j++) { 912 row = *adj++; 913 pnz = pi_loc[row+1] - pi_loc[row]; 914 Jptr = pj_loc + pi_loc[row]; 915 /* Expand list if it is not long enough */ 916 if (pnz+apnz_max > lsize) { 917 lsize = pnz+apnz_max; 918 ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 919 } 920 /* add non-zero cols of P into the sorted linked list lnk */ 921 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 922 apnz = *lnk; /* The first element in the list is the number of items in the list */ 923 api[i+1] = api[i] + apnz; 924 if (apnz > apnz_max) apnz_max = apnz; 925 } 926 /* off-diagonal portion of A */ 927 nzi = aoi[i+1] - aoi[i]; 928 for (j=0; j<nzi; j++) { 929 row = *aoj++; 930 pnz = pi_oth[row+1] - pi_oth[row]; 931 Jptr = pj_oth + pi_oth[row]; 932 /* Expand list if it is not long enough */ 933 if (pnz+apnz_max > lsize) { 934 lsize = pnz + apnz_max; 935 ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 936 } 937 /* add non-zero cols of P into the sorted linked list lnk */ 938 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 939 apnz = *lnk; /* The first element in the list is the number of items in the list */ 940 api[i+1] = api[i] + apnz; 941 if (apnz > apnz_max) apnz_max = apnz; 942 } 943 apnz = *lnk; 944 api[i+1] = api[i] + apnz; 945 if (apnz > apnz_max) apnz_max = apnz; 946 947 /* if free space is not available, double the total space in the list */ 948 if (current_space->local_remaining<apnz) { 949 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 950 nspacedouble++; 951 } 952 953 /* Copy data into free space, then initialize lnk */ 954 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 955 ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 956 957 current_space->array += apnz; 958 current_space->local_used += apnz; 959 current_space->local_remaining -= apnz; 960 } 961 962 /* Allocate space for apj, initialize apj, and */ 963 /* destroy list of free space and other temporary array(s) */ 964 ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 965 apj = ptap->apj; 966 ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 967 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 968 969 /* create and assemble symbolic parallel matrix Cmpi */ 970 /*----------------------------------------------------*/ 971 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 972 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 973 ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 974 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 975 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 976 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 977 978 /* malloc apa for assembly Cmpi */ 979 ierr = PetscCalloc1(apnz_max,&ptap->apa);CHKERRQ(ierr); 980 981 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 982 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 983 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 984 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 985 986 ptap->destroy = Cmpi->ops->destroy; 987 ptap->duplicate = Cmpi->ops->duplicate; 988 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 989 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 990 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 991 992 /* attach the supporting struct to Cmpi for reuse */ 993 c = (Mat_MPIAIJ*)Cmpi->data; 994 c->ap = ptap; 995 *C = Cmpi; 996 997 /* set MatInfo */ 998 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 999 if (afill < 1.0) afill = 1.0; 1000 Cmpi->info.mallocs = nspacedouble; 1001 Cmpi->info.fill_ratio_given = fill; 1002 Cmpi->info.fill_ratio_needed = afill; 1003 1004 #if defined(PETSC_USE_INFO) 1005 if (api[am]) { 1006 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1007 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 1008 } else { 1009 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1010 } 1011 #endif 1012 PetscFunctionReturn(0); 1013 } 1014 1015 /* This function is needed for the seqMPI matrix-matrix multiplication. */ 1016 /* Three input arrays are merged to one output array. The size of the */ 1017 /* output array is also output. Duplicate entries only show up once. */ 1018 static void Merge3SortedArrays(PetscInt size1, PetscInt *in1, 1019 PetscInt size2, PetscInt *in2, 1020 PetscInt size3, PetscInt *in3, 1021 PetscInt *size4, PetscInt *out) 1022 { 1023 int i = 0, j = 0, k = 0, l = 0; 1024 1025 /* Traverse all three arrays */ 1026 while (i<size1 && j<size2 && k<size3) { 1027 if (in1[i] < in2[j] && in1[i] < in3[k]) { 1028 out[l++] = in1[i++]; 1029 } 1030 else if(in2[j] < in1[i] && in2[j] < in3[k]) { 1031 out[l++] = in2[j++]; 1032 } 1033 else if(in3[k] < in1[i] && in3[k] < in2[j]) { 1034 out[l++] = in3[k++]; 1035 } 1036 else if(in1[i] == in2[j] && in1[i] < in3[k]) { 1037 out[l++] = in1[i]; 1038 i++, j++; 1039 } 1040 else if(in1[i] == in3[k] && in1[i] < in2[j]) { 1041 out[l++] = in1[i]; 1042 i++, k++; 1043 } 1044 else if(in3[k] == in2[j] && in2[j] < in1[i]) { 1045 out[l++] = in2[j]; 1046 k++, j++; 1047 } 1048 else if(in1[i] == in2[j] && in1[i] == in3[k]) { 1049 out[l++] = in1[i]; 1050 i++, j++, k++; 1051 } 1052 } 1053 1054 /* Traverse two remaining arrays */ 1055 while (i<size1 && j<size2) { 1056 if (in1[i] < in2[j]) { 1057 out[l++] = in1[i++]; 1058 } 1059 else if(in1[i] > in2[j]) { 1060 out[l++] = in2[j++]; 1061 } 1062 else { 1063 out[l++] = in1[i]; 1064 i++, j++; 1065 } 1066 } 1067 1068 while (i<size1 && k<size3) { 1069 if (in1[i] < in3[k]) { 1070 out[l++] = in1[i++]; 1071 } 1072 else if(in1[i] > in3[k]) { 1073 out[l++] = in3[k++]; 1074 } 1075 else { 1076 out[l++] = in1[i]; 1077 i++, k++; 1078 } 1079 } 1080 1081 while (k<size3 && j<size2) { 1082 if (in3[k] < in2[j]) { 1083 out[l++] = in3[k++]; 1084 } 1085 else if(in3[k] > in2[j]) { 1086 out[l++] = in2[j++]; 1087 } 1088 else { 1089 out[l++] = in3[k]; 1090 k++, j++; 1091 } 1092 } 1093 1094 /* Traverse one remaining array */ 1095 while (i<size1) out[l++] = in1[i++]; 1096 while (j<size2) out[l++] = in2[j++]; 1097 while (k<size3) out[l++] = in3[k++]; 1098 1099 *size4 = l; 1100 } 1101 1102 /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and */ 1103 /* adds up the products. Two of these three multiplications are performed with existing (sequential) */ 1104 /* matrix-matrix multiplications. */ 1105 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat *C) 1106 { 1107 PetscErrorCode ierr; 1108 MPI_Comm comm; 1109 PetscMPIInt size; 1110 Mat Cmpi; 1111 Mat_APMPI *ptap; 1112 PetscFreeSpaceList free_space_diag=NULL, current_space=NULL; 1113 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data; 1114 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc; 1115 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data; 1116 Mat_MPIAIJ *c; 1117 Mat_SeqAIJ *adpd_seq, *p_off, *aopoth_seq; 1118 PetscInt adponz, adpdnz; 1119 PetscInt *pi_loc,*dnz,*onz; 1120 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart; 1121 PetscInt *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi, 1122 *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj; 1123 PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend; 1124 PetscBT lnkbt; 1125 PetscReal afill; 1126 PetscMPIInt rank; 1127 Mat adpd, aopoth; 1128 MatType mtype; 1129 const char *prefix; 1130 1131 PetscFunctionBegin; 1132 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1133 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1134 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1135 ierr = MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend); CHKERRQ(ierr); 1136 1137 /* create struct Mat_APMPI and attached it to C later */ 1138 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1139 1140 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1141 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1142 1143 /* get P_loc by taking all local rows of P */ 1144 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1145 1146 1147 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1148 pi_loc = p_loc->i; 1149 1150 /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */ 1151 ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 1152 ierr = PetscMalloc1(am+2,&adpoi);CHKERRQ(ierr); 1153 1154 adpoi[0] = 0; 1155 ptap->api = api; 1156 api[0] = 0; 1157 1158 /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */ 1159 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1160 ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 1161 1162 /* Symbolic calc of A_loc_diag * P_loc_diag */ 1163 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1164 ierr = MatSetOptionsPrefix(a->A,prefix);CHKERRQ(ierr); 1165 ierr = MatAppendOptionsPrefix(a->A,"inner_diag_");CHKERRQ(ierr); 1166 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);CHKERRQ(ierr); 1167 adpd_seq = (Mat_SeqAIJ*)((adpd)->data); 1168 adpdi = adpd_seq->i; adpdj = adpd_seq->j; 1169 p_off = (Mat_SeqAIJ*)((p->B)->data); 1170 poff_i = p_off->i; poff_j = p_off->j; 1171 1172 /* j_temp stores indices of a result row before they are added to the linked list */ 1173 ierr = PetscMalloc1(pN+2,&j_temp);CHKERRQ(ierr); 1174 1175 1176 /* Symbolic calc of the A_diag * p_loc_off */ 1177 /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 1178 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);CHKERRQ(ierr); 1179 current_space = free_space_diag; 1180 1181 for (i=0; i<am; i++) { 1182 /* A_diag * P_loc_off */ 1183 nzi = adi[i+1] - adi[i]; 1184 for (j=0; j<nzi; j++) { 1185 row = *adj++; 1186 pnz = poff_i[row+1] - poff_i[row]; 1187 Jptr = poff_j + poff_i[row]; 1188 for(i1 = 0; i1 < pnz; i1++) { 1189 j_temp[i1] = p->garray[Jptr[i1]]; 1190 } 1191 /* add non-zero cols of P into the sorted linked list lnk */ 1192 ierr = PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);CHKERRQ(ierr); 1193 } 1194 1195 adponz = lnk[0]; 1196 adpoi[i+1] = adpoi[i] + adponz; 1197 1198 /* if free space is not available, double the total space in the list */ 1199 if (current_space->local_remaining<adponz) { 1200 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1201 nspacedouble++; 1202 } 1203 1204 /* Copy data into free space, then initialize lnk */ 1205 ierr = PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1206 1207 current_space->array += adponz; 1208 current_space->local_used += adponz; 1209 current_space->local_remaining -= adponz; 1210 } 1211 1212 /* Symbolic calc of A_off * P_oth */ 1213 ierr = MatSetOptionsPrefix(a->B,prefix);CHKERRQ(ierr); 1214 ierr = MatAppendOptionsPrefix(a->B,"inner_offdiag_");CHKERRQ(ierr); 1215 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);CHKERRQ(ierr); 1216 aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data); 1217 aopothi = aopoth_seq->i; aopothj = aopoth_seq->j; 1218 1219 /* Allocate space for apj, adpj, aopj, ... */ 1220 /* destroy lists of free space and other temporary array(s) */ 1221 1222 ierr = PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);CHKERRQ(ierr); 1223 ierr = PetscMalloc1(adpoi[am]+2, &adpoj);CHKERRQ(ierr); 1224 1225 /* Copy from linked list to j-array */ 1226 ierr = PetscFreeSpaceContiguous(&free_space_diag,adpoj);CHKERRQ(ierr); 1227 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1228 1229 adpoJ = adpoj; 1230 adpdJ = adpdj; 1231 aopJ = aopothj; 1232 apj = ptap->apj; 1233 apJ = apj; /* still empty */ 1234 1235 /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */ 1236 /* A_diag * P_loc_diag to get A*P */ 1237 for (i = 0; i < am; i++) { 1238 aopnz = aopothi[i+1] - aopothi[i]; 1239 adponz = adpoi[i+1] - adpoi[i]; 1240 adpdnz = adpdi[i+1] - adpdi[i]; 1241 1242 /* Correct indices from A_diag*P_diag */ 1243 for(i1 = 0; i1 < adpdnz; i1++) { 1244 adpdJ[i1] += p_colstart; 1245 } 1246 /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */ 1247 Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ); 1248 ierr = MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz); CHKERRQ(ierr); 1249 1250 aopJ += aopnz; 1251 adpoJ += adponz; 1252 adpdJ += adpdnz; 1253 apJ += apnz; 1254 api[i+1] = api[i] + apnz; 1255 } 1256 1257 /* malloc apa to store dense row A[i,:]*P */ 1258 ierr = PetscCalloc1(pN+2,&ptap->apa);CHKERRQ(ierr); 1259 1260 /* create and assemble symbolic parallel matrix Cmpi */ 1261 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1262 ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1263 ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 1264 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1265 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1266 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1267 1268 1269 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 1270 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1271 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1272 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1273 1274 1275 ptap->destroy = Cmpi->ops->destroy; 1276 ptap->duplicate = Cmpi->ops->duplicate; 1277 Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1278 Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1279 1280 /* attach the supporting struct to Cmpi for reuse */ 1281 c = (Mat_MPIAIJ*)Cmpi->data; 1282 c->ap = ptap; 1283 *C = Cmpi; 1284 1285 /* set MatInfo */ 1286 afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 1287 if (afill < 1.0) afill = 1.0; 1288 Cmpi->info.mallocs = nspacedouble; 1289 Cmpi->info.fill_ratio_given = fill; 1290 Cmpi->info.fill_ratio_needed = afill; 1291 1292 #if defined(PETSC_USE_INFO) 1293 if (api[am]) { 1294 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1295 ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 1296 } else { 1297 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1298 } 1299 #endif 1300 1301 ierr = MatDestroy(&aopoth);CHKERRQ(ierr); 1302 ierr = MatDestroy(&adpd);CHKERRQ(ierr); 1303 ierr = PetscFree(j_temp);CHKERRQ(ierr); 1304 ierr = PetscFree(adpoj);CHKERRQ(ierr); 1305 ierr = PetscFree(adpoi);CHKERRQ(ierr); 1306 PetscFunctionReturn(0); 1307 } 1308 1309 1310 /*-------------------------------------------------------------------------*/ 1311 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 1312 { 1313 PetscErrorCode ierr; 1314 const char *algTypes[3] = {"scalable","nonscalable","matmatmult"}; 1315 PetscInt aN=A->cmap->N,alg=1; /* set default algorithm */ 1316 PetscBool flg; 1317 1318 PetscFunctionBegin; 1319 if (scall == MAT_INITIAL_MATRIX) { 1320 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr); 1321 ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[1],&alg,&flg);CHKERRQ(ierr); 1322 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1323 1324 ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 1325 switch (alg) { 1326 case 1: 1327 if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */ 1328 MatInfo Ainfo,Pinfo; 1329 PetscInt nz_local; 1330 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 1331 MPI_Comm comm; 1332 1333 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 1334 ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 1335 nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); /* estimated local nonzero entries */ 1336 1337 if (aN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 1338 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1339 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 1340 1341 if (alg_scalable) { 1342 alg = 0; /* scalable algorithm would slower than nonscalable algorithm */ 1343 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1344 break; 1345 } 1346 } 1347 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr); 1348 break; 1349 case 2: 1350 { 1351 Mat Pt; 1352 Mat_APMPI *ptap; 1353 Mat_MPIAIJ *c; 1354 ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr); 1355 ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr); 1356 c = (Mat_MPIAIJ*)(*C)->data; 1357 ptap = c->ap; 1358 if (ptap) { 1359 ptap->Pt = Pt; 1360 (*C)->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1361 } 1362 (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult; 1363 PetscFunctionReturn(0); 1364 } 1365 break; 1366 default: /* scalable algorithm */ 1367 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1368 break; 1369 } 1370 ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 1371 1372 { 1373 Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 1374 Mat_APMPI *ap = c->ap; 1375 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 1376 ap->freestruct = PETSC_FALSE; 1377 ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 1378 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1379 } 1380 } 1381 1382 ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 1383 ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 1384 ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 1385 PetscFunctionReturn(0); 1386 } 1387 1388 /* This routine only works when scall=MAT_REUSE_MATRIX! */ 1389 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C) 1390 { 1391 PetscErrorCode ierr; 1392 Mat_MPIAIJ *c=(Mat_MPIAIJ*)C->data; 1393 Mat_APMPI *ptap= c->ap; 1394 Mat Pt; 1395 1396 PetscFunctionBegin; 1397 if (!ptap->Pt) { 1398 MPI_Comm comm; 1399 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1400 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1401 } 1402 1403 Pt=ptap->Pt; 1404 ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr); 1405 ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr); 1406 1407 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1408 if (ptap->freestruct) { 1409 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 1410 } 1411 PetscFunctionReturn(0); 1412 } 1413 1414 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1415 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C) 1416 { 1417 PetscErrorCode ierr; 1418 Mat_APMPI *ptap; 1419 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1420 MPI_Comm comm; 1421 PetscMPIInt size,rank; 1422 Mat Cmpi; 1423 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1424 PetscInt pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n; 1425 PetscInt *lnk,i,k,nsend; 1426 PetscBT lnkbt; 1427 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 1428 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1429 PetscInt len,proc,*dnz,*onz,*owners,nzi; 1430 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1431 MPI_Request *swaits,*rwaits; 1432 MPI_Status *sstatus,rstatus; 1433 PetscLayout rowmap; 1434 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 1435 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 1436 PetscInt *Jptr,*prmap=p->garray,con,j,Crmax; 1437 Mat_SeqAIJ *a_loc,*c_loc,*c_oth; 1438 PetscTable ta; 1439 MatType mtype; 1440 const char *prefix; 1441 1442 PetscFunctionBegin; 1443 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1444 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1445 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1446 1447 /* create symbolic parallel matrix Cmpi */ 1448 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1449 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1450 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1451 1452 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1453 1454 /* create struct Mat_APMPI and attached it to C later */ 1455 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1456 ptap->reuse = MAT_INITIAL_MATRIX; 1457 1458 /* (0) compute Rd = Pd^T, Ro = Po^T */ 1459 /* --------------------------------- */ 1460 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1461 ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1462 1463 /* (1) compute symbolic A_loc */ 1464 /* ---------------------------*/ 1465 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1466 1467 /* (2-1) compute symbolic C_oth = Ro*A_loc */ 1468 /* ------------------------------------ */ 1469 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1470 ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 1471 ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 1472 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 1473 1474 /* (3) send coj of C_oth to other processors */ 1475 /* ------------------------------------------ */ 1476 /* determine row ownership */ 1477 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 1478 rowmap->n = pn; 1479 rowmap->bs = 1; 1480 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 1481 owners = rowmap->range; 1482 1483 /* determine the number of messages to send, their lengths */ 1484 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 1485 ierr = PetscArrayzero(len_s,size);CHKERRQ(ierr); 1486 ierr = PetscArrayzero(len_si,size);CHKERRQ(ierr); 1487 1488 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 1489 coi = c_oth->i; coj = c_oth->j; 1490 con = ptap->C_oth->rmap->n; 1491 proc = 0; 1492 for (i=0; i<con; i++) { 1493 while (prmap[i] >= owners[proc+1]) proc++; 1494 len_si[proc]++; /* num of rows in Co(=Pt*A) to be sent to [proc] */ 1495 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 1496 } 1497 1498 len = 0; /* max length of buf_si[], see (4) */ 1499 owners_co[0] = 0; 1500 nsend = 0; 1501 for (proc=0; proc<size; proc++) { 1502 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1503 if (len_s[proc]) { 1504 nsend++; 1505 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 1506 len += len_si[proc]; 1507 } 1508 } 1509 1510 /* determine the number and length of messages to receive for coi and coj */ 1511 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 1512 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 1513 1514 /* post the Irecv and Isend of coj */ 1515 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1516 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1517 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 1518 for (proc=0, k=0; proc<size; proc++) { 1519 if (!len_s[proc]) continue; 1520 i = owners_co[proc]; 1521 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1522 k++; 1523 } 1524 1525 /* (2-2) compute symbolic C_loc = Rd*A_loc */ 1526 /* ---------------------------------------- */ 1527 ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 1528 ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 1529 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 1530 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 1531 1532 /* receives coj are complete */ 1533 for (i=0; i<nrecv; i++) { 1534 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1535 } 1536 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1537 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1538 1539 /* add received column indices into ta to update Crmax */ 1540 a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data; 1541 1542 /* create and initialize a linked list */ 1543 ierr = PetscTableCreate(an,aN,&ta);CHKERRQ(ierr); /* for compute Crmax */ 1544 MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta); 1545 1546 for (k=0; k<nrecv; k++) {/* k-th received message */ 1547 Jptr = buf_rj[k]; 1548 for (j=0; j<len_r[k]; j++) { 1549 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1550 } 1551 } 1552 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 1553 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 1554 1555 /* (4) send and recv coi */ 1556 /*-----------------------*/ 1557 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1558 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1559 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1560 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1561 for (proc=0,k=0; proc<size; proc++) { 1562 if (!len_s[proc]) continue; 1563 /* form outgoing message for i-structure: 1564 buf_si[0]: nrows to be sent 1565 [1:nrows]: row index (global) 1566 [nrows+1:2*nrows+1]: i-structure index 1567 */ 1568 /*-------------------------------------------*/ 1569 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 1570 buf_si_i = buf_si + nrows+1; 1571 buf_si[0] = nrows; 1572 buf_si_i[0] = 0; 1573 nrows = 0; 1574 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1575 nzi = coi[i+1] - coi[i]; 1576 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1577 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1578 nrows++; 1579 } 1580 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1581 k++; 1582 buf_si += len_si[proc]; 1583 } 1584 for (i=0; i<nrecv; i++) { 1585 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1586 } 1587 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1588 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1589 1590 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 1591 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1592 ierr = PetscFree(swaits);CHKERRQ(ierr); 1593 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1594 1595 /* (5) compute the local portion of Cmpi */ 1596 /* ------------------------------------------ */ 1597 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 1598 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 1599 current_space = free_space; 1600 1601 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 1602 for (k=0; k<nrecv; k++) { 1603 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1604 nrows = *buf_ri_k[k]; 1605 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1606 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1607 } 1608 1609 ierr = MatPreallocateInitialize(comm,pn,an,dnz,onz);CHKERRQ(ierr); 1610 ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1611 for (i=0; i<pn; i++) { 1612 /* add C_loc into Cmpi */ 1613 nzi = c_loc->i[i+1] - c_loc->i[i]; 1614 Jptr = c_loc->j + c_loc->i[i]; 1615 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1616 1617 /* add received col data into lnk */ 1618 for (k=0; k<nrecv; k++) { /* k-th received message */ 1619 if (i == *nextrow[k]) { /* i-th row */ 1620 nzi = *(nextci[k]+1) - *nextci[k]; 1621 Jptr = buf_rj[k] + *nextci[k]; 1622 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1623 nextrow[k]++; nextci[k]++; 1624 } 1625 } 1626 nzi = lnk[0]; 1627 1628 /* copy data into free space, then initialize lnk */ 1629 ierr = PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1630 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 1631 } 1632 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1633 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1634 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 1635 1636 /* local sizes and preallocation */ 1637 ierr = MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1638 if (P->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr);} 1639 if (A->cmap->bs > 0) {ierr = PetscLayoutSetBlockSize(Cmpi->cmap,A->cmap->bs);CHKERRQ(ierr);} 1640 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1641 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1642 1643 /* members in merge */ 1644 ierr = PetscFree(id_r);CHKERRQ(ierr); 1645 ierr = PetscFree(len_r);CHKERRQ(ierr); 1646 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 1647 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 1648 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 1649 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 1650 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 1651 1652 /* attach the supporting struct to Cmpi for reuse */ 1653 c = (Mat_MPIAIJ*)Cmpi->data; 1654 c->ap = ptap; 1655 ptap->destroy = Cmpi->ops->destroy; 1656 1657 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1658 Cmpi->assembled = PETSC_FALSE; 1659 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1660 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1661 1662 *C = Cmpi; 1663 PetscFunctionReturn(0); 1664 } 1665 1666 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C) 1667 { 1668 PetscErrorCode ierr; 1669 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1670 Mat_SeqAIJ *c_seq; 1671 Mat_APMPI *ptap = c->ap; 1672 Mat A_loc,C_loc,C_oth; 1673 PetscInt i,rstart,rend,cm,ncols,row; 1674 const PetscInt *cols; 1675 const PetscScalar *vals; 1676 1677 PetscFunctionBegin; 1678 if (!ptap->A_loc) { 1679 MPI_Comm comm; 1680 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1681 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1682 } 1683 1684 ierr = MatZeroEntries(C);CHKERRQ(ierr); 1685 1686 if (ptap->reuse == MAT_REUSE_MATRIX) { 1687 /* These matrices are obtained in MatTransposeMatMultSymbolic() */ 1688 /* 1) get R = Pd^T, Ro = Po^T */ 1689 /*----------------------------*/ 1690 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1691 ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1692 1693 /* 2) compute numeric A_loc */ 1694 /*--------------------------*/ 1695 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1696 } 1697 1698 /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */ 1699 A_loc = ptap->A_loc; 1700 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);CHKERRQ(ierr); 1701 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);CHKERRQ(ierr); 1702 C_loc = ptap->C_loc; 1703 C_oth = ptap->C_oth; 1704 1705 /* add C_loc and Co to to C */ 1706 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 1707 1708 /* C_loc -> C */ 1709 cm = C_loc->rmap->N; 1710 c_seq = (Mat_SeqAIJ*)C_loc->data; 1711 cols = c_seq->j; 1712 vals = c_seq->a; 1713 for (i=0; i<cm; i++) { 1714 ncols = c_seq->i[i+1] - c_seq->i[i]; 1715 row = rstart + i; 1716 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1717 cols += ncols; vals += ncols; 1718 } 1719 1720 /* Co -> C, off-processor part */ 1721 cm = C_oth->rmap->N; 1722 c_seq = (Mat_SeqAIJ*)C_oth->data; 1723 cols = c_seq->j; 1724 vals = c_seq->a; 1725 for (i=0; i<cm; i++) { 1726 ncols = c_seq->i[i+1] - c_seq->i[i]; 1727 row = p->garray[i]; 1728 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1729 cols += ncols; vals += ncols; 1730 } 1731 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1732 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1733 1734 ptap->reuse = MAT_REUSE_MATRIX; 1735 1736 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1737 if (ptap->freestruct) { 1738 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 1739 } 1740 PetscFunctionReturn(0); 1741 } 1742 1743 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1744 { 1745 PetscErrorCode ierr; 1746 Mat_Merge_SeqsToMPI *merge; 1747 Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1748 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1749 Mat_APMPI *ptap; 1750 PetscInt *adj; 1751 PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1752 MatScalar *ada,*ca,valtmp; 1753 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1754 MPI_Comm comm; 1755 PetscMPIInt size,rank,taga,*len_s; 1756 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1757 PetscInt **buf_ri,**buf_rj; 1758 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1759 MPI_Request *s_waits,*r_waits; 1760 MPI_Status *status; 1761 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1762 PetscInt *ai,*aj,*coi,*coj,*poJ,*pdJ; 1763 Mat A_loc; 1764 Mat_SeqAIJ *a_loc; 1765 1766 PetscFunctionBegin; 1767 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1768 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1769 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1770 1771 ptap = c->ap; 1772 if (!ptap->A_loc) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1773 merge = ptap->merge; 1774 1775 /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1776 /*------------------------------------------*/ 1777 /* get data from symbolic products */ 1778 coi = merge->coi; coj = merge->coj; 1779 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1780 bi = merge->bi; bj = merge->bj; 1781 owners = merge->rowmap->range; 1782 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1783 1784 /* get A_loc by taking all local rows of A */ 1785 A_loc = ptap->A_loc; 1786 ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1787 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1788 ai = a_loc->i; 1789 aj = a_loc->j; 1790 1791 for (i=0; i<am; i++) { 1792 anz = ai[i+1] - ai[i]; 1793 adj = aj + ai[i]; 1794 ada = a_loc->a + ai[i]; 1795 1796 /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1797 /*-------------------------------------------------------------*/ 1798 /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1799 pnz = po->i[i+1] - po->i[i]; 1800 poJ = po->j + po->i[i]; 1801 pA = po->a + po->i[i]; 1802 for (j=0; j<pnz; j++) { 1803 row = poJ[j]; 1804 cj = coj + coi[row]; 1805 ca = coa + coi[row]; 1806 /* perform sparse axpy */ 1807 nexta = 0; 1808 valtmp = pA[j]; 1809 for (k=0; nexta<anz; k++) { 1810 if (cj[k] == adj[nexta]) { 1811 ca[k] += valtmp*ada[nexta]; 1812 nexta++; 1813 } 1814 } 1815 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1816 } 1817 1818 /* put the value into Cd (diagonal part) */ 1819 pnz = pd->i[i+1] - pd->i[i]; 1820 pdJ = pd->j + pd->i[i]; 1821 pA = pd->a + pd->i[i]; 1822 for (j=0; j<pnz; j++) { 1823 row = pdJ[j]; 1824 cj = bj + bi[row]; 1825 ca = ba + bi[row]; 1826 /* perform sparse axpy */ 1827 nexta = 0; 1828 valtmp = pA[j]; 1829 for (k=0; nexta<anz; k++) { 1830 if (cj[k] == adj[nexta]) { 1831 ca[k] += valtmp*ada[nexta]; 1832 nexta++; 1833 } 1834 } 1835 ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1836 } 1837 } 1838 1839 /* 3) send and recv matrix values coa */ 1840 /*------------------------------------*/ 1841 buf_ri = merge->buf_ri; 1842 buf_rj = merge->buf_rj; 1843 len_s = merge->len_s; 1844 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1845 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1846 1847 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1848 for (proc=0,k=0; proc<size; proc++) { 1849 if (!len_s[proc]) continue; 1850 i = merge->owners_co[proc]; 1851 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1852 k++; 1853 } 1854 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1855 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1856 1857 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1858 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1859 ierr = PetscFree(coa);CHKERRQ(ierr); 1860 1861 /* 4) insert local Cseq and received values into Cmpi */ 1862 /*----------------------------------------------------*/ 1863 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1864 for (k=0; k<merge->nrecv; k++) { 1865 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1866 nrows = *(buf_ri_k[k]); 1867 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1868 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1869 } 1870 1871 for (i=0; i<cm; i++) { 1872 row = owners[rank] + i; /* global row index of C_seq */ 1873 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1874 ba_i = ba + bi[i]; 1875 bnz = bi[i+1] - bi[i]; 1876 /* add received vals into ba */ 1877 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1878 /* i-th row */ 1879 if (i == *nextrow[k]) { 1880 cnz = *(nextci[k]+1) - *nextci[k]; 1881 cj = buf_rj[k] + *(nextci[k]); 1882 ca = abuf_r[k] + *(nextci[k]); 1883 nextcj = 0; 1884 for (j=0; nextcj<cnz; j++) { 1885 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1886 ba_i[j] += ca[nextcj++]; 1887 } 1888 } 1889 nextrow[k]++; nextci[k]++; 1890 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1891 } 1892 } 1893 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1894 } 1895 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1896 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1897 1898 ierr = PetscFree(ba);CHKERRQ(ierr); 1899 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1900 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1901 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1902 1903 if (ptap->freestruct) { 1904 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 1905 } 1906 PetscFunctionReturn(0); 1907 } 1908 1909 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1910 { 1911 PetscErrorCode ierr; 1912 Mat Cmpi,A_loc,POt,PDt; 1913 Mat_APMPI *ptap; 1914 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1915 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c; 1916 PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1917 PetscInt nnz; 1918 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1919 PetscInt am =A->rmap->n,pn=P->cmap->n; 1920 MPI_Comm comm; 1921 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1922 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1923 PetscInt len,proc,*dnz,*onz,*owners; 1924 PetscInt nzi,*bi,*bj; 1925 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1926 MPI_Request *swaits,*rwaits; 1927 MPI_Status *sstatus,rstatus; 1928 Mat_Merge_SeqsToMPI *merge; 1929 PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1930 PetscReal afill =1.0,afill_tmp; 1931 PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax; 1932 Mat_SeqAIJ *a_loc,*pdt,*pot; 1933 PetscTable ta; 1934 MatType mtype; 1935 1936 PetscFunctionBegin; 1937 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1938 /* check if matrix local sizes are compatible */ 1939 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); 1940 1941 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1942 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1943 1944 /* create struct Mat_APMPI and attached it to C later */ 1945 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1946 1947 /* get A_loc by taking all local rows of A */ 1948 ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1949 1950 ptap->A_loc = A_loc; 1951 a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1952 ai = a_loc->i; 1953 aj = a_loc->j; 1954 1955 /* determine symbolic Co=(p->B)^T*A - send to others */ 1956 /*----------------------------------------------------*/ 1957 ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1958 pdt = (Mat_SeqAIJ*)PDt->data; 1959 pdti = pdt->i; pdtj = pdt->j; 1960 1961 ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1962 pot = (Mat_SeqAIJ*)POt->data; 1963 poti = pot->i; potj = pot->j; 1964 1965 /* then, compute symbolic Co = (p->B)^T*A */ 1966 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1967 >= (num of nonzero rows of C_seq) - pn */ 1968 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1969 coi[0] = 0; 1970 1971 /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1972 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am])); 1973 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1974 current_space = free_space; 1975 1976 /* create and initialize a linked list */ 1977 ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr); 1978 MatRowMergeMax_SeqAIJ(a_loc,am,ta); 1979 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 1980 1981 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 1982 1983 for (i=0; i<pon; i++) { 1984 pnz = poti[i+1] - poti[i]; 1985 ptJ = potj + poti[i]; 1986 for (j=0; j<pnz; j++) { 1987 row = ptJ[j]; /* row of A_loc == col of Pot */ 1988 anz = ai[row+1] - ai[row]; 1989 Jptr = aj + ai[row]; 1990 /* add non-zero cols of AP into the sorted linked list lnk */ 1991 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1992 } 1993 nnz = lnk[0]; 1994 1995 /* If free space is not available, double the total space in the list */ 1996 if (current_space->local_remaining<nnz) { 1997 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1998 nspacedouble++; 1999 } 2000 2001 /* Copy data into free space, and zero out denserows */ 2002 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2003 2004 current_space->array += nnz; 2005 current_space->local_used += nnz; 2006 current_space->local_remaining -= nnz; 2007 2008 coi[i+1] = coi[i] + nnz; 2009 } 2010 2011 ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 2012 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 2013 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */ 2014 2015 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 2016 if (afill_tmp > afill) afill = afill_tmp; 2017 2018 /* send j-array (coj) of Co to other processors */ 2019 /*----------------------------------------------*/ 2020 /* determine row ownership */ 2021 ierr = PetscNew(&merge);CHKERRQ(ierr); 2022 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 2023 2024 merge->rowmap->n = pn; 2025 merge->rowmap->bs = 1; 2026 2027 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 2028 owners = merge->rowmap->range; 2029 2030 /* determine the number of messages to send, their lengths */ 2031 ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 2032 ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 2033 2034 len_s = merge->len_s; 2035 merge->nsend = 0; 2036 2037 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 2038 2039 proc = 0; 2040 for (i=0; i<pon; i++) { 2041 while (prmap[i] >= owners[proc+1]) proc++; 2042 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 2043 len_s[proc] += coi[i+1] - coi[i]; 2044 } 2045 2046 len = 0; /* max length of buf_si[] */ 2047 owners_co[0] = 0; 2048 for (proc=0; proc<size; proc++) { 2049 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 2050 if (len_si[proc]) { 2051 merge->nsend++; 2052 len_si[proc] = 2*(len_si[proc] + 1); 2053 len += len_si[proc]; 2054 } 2055 } 2056 2057 /* determine the number and length of messages to receive for coi and coj */ 2058 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 2059 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 2060 2061 /* post the Irecv and Isend of coj */ 2062 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 2063 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 2064 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 2065 for (proc=0, k=0; proc<size; proc++) { 2066 if (!len_s[proc]) continue; 2067 i = owners_co[proc]; 2068 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 2069 k++; 2070 } 2071 2072 /* receives and sends of coj are complete */ 2073 ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 2074 for (i=0; i<merge->nrecv; i++) { 2075 PetscMPIInt icompleted; 2076 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 2077 } 2078 ierr = PetscFree(rwaits);CHKERRQ(ierr); 2079 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 2080 2081 /* add received column indices into table to update Armax */ 2082 /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */ 2083 for (k=0; k<merge->nrecv; k++) {/* k-th received message */ 2084 Jptr = buf_rj[k]; 2085 for (j=0; j<merge->len_r[k]; j++) { 2086 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 2087 } 2088 } 2089 ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 2090 /* printf("Armax %d, an %d + Bn %d = %d, aN %d\n",Armax,A->cmap->n,a->B->cmap->N,A->cmap->n+a->B->cmap->N,aN); */ 2091 2092 /* send and recv coi */ 2093 /*-------------------*/ 2094 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 2095 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 2096 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 2097 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 2098 for (proc=0,k=0; proc<size; proc++) { 2099 if (!len_s[proc]) continue; 2100 /* form outgoing message for i-structure: 2101 buf_si[0]: nrows to be sent 2102 [1:nrows]: row index (global) 2103 [nrows+1:2*nrows+1]: i-structure index 2104 */ 2105 /*-------------------------------------------*/ 2106 nrows = len_si[proc]/2 - 1; 2107 buf_si_i = buf_si + nrows+1; 2108 buf_si[0] = nrows; 2109 buf_si_i[0] = 0; 2110 nrows = 0; 2111 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 2112 nzi = coi[i+1] - coi[i]; 2113 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 2114 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 2115 nrows++; 2116 } 2117 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 2118 k++; 2119 buf_si += len_si[proc]; 2120 } 2121 i = merge->nrecv; 2122 while (i--) { 2123 PetscMPIInt icompleted; 2124 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 2125 } 2126 ierr = PetscFree(rwaits);CHKERRQ(ierr); 2127 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 2128 ierr = PetscFree(len_si);CHKERRQ(ierr); 2129 ierr = PetscFree(len_ri);CHKERRQ(ierr); 2130 ierr = PetscFree(swaits);CHKERRQ(ierr); 2131 ierr = PetscFree(sstatus);CHKERRQ(ierr); 2132 ierr = PetscFree(buf_s);CHKERRQ(ierr); 2133 2134 /* compute the local portion of C (mpi mat) */ 2135 /*------------------------------------------*/ 2136 /* allocate bi array and free space for accumulating nonzero column info */ 2137 ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 2138 bi[0] = 0; 2139 2140 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 2141 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am]))); 2142 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 2143 current_space = free_space; 2144 2145 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 2146 for (k=0; k<merge->nrecv; k++) { 2147 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 2148 nrows = *buf_ri_k[k]; 2149 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 2150 nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure */ 2151 } 2152 2153 ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 2154 ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 2155 rmax = 0; 2156 for (i=0; i<pn; i++) { 2157 /* add pdt[i,:]*AP into lnk */ 2158 pnz = pdti[i+1] - pdti[i]; 2159 ptJ = pdtj + pdti[i]; 2160 for (j=0; j<pnz; j++) { 2161 row = ptJ[j]; /* row of AP == col of Pt */ 2162 anz = ai[row+1] - ai[row]; 2163 Jptr = aj + ai[row]; 2164 /* add non-zero cols of AP into the sorted linked list lnk */ 2165 ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 2166 } 2167 2168 /* add received col data into lnk */ 2169 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 2170 if (i == *nextrow[k]) { /* i-th row */ 2171 nzi = *(nextci[k]+1) - *nextci[k]; 2172 Jptr = buf_rj[k] + *nextci[k]; 2173 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 2174 nextrow[k]++; nextci[k]++; 2175 } 2176 } 2177 nnz = lnk[0]; 2178 2179 /* if free space is not available, make more free space */ 2180 if (current_space->local_remaining<nnz) { 2181 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 2182 nspacedouble++; 2183 } 2184 /* copy data into free space, then initialize lnk */ 2185 ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2186 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2187 2188 current_space->array += nnz; 2189 current_space->local_used += nnz; 2190 current_space->local_remaining -= nnz; 2191 2192 bi[i+1] = bi[i] + nnz; 2193 if (nnz > rmax) rmax = nnz; 2194 } 2195 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 2196 2197 ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 2198 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2199 afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 2200 if (afill_tmp > afill) afill = afill_tmp; 2201 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 2202 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 2203 2204 ierr = MatDestroy(&POt);CHKERRQ(ierr); 2205 ierr = MatDestroy(&PDt);CHKERRQ(ierr); 2206 2207 /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 2208 /*----------------------------------------------------------------------------------*/ 2209 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 2210 ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 2211 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 2212 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 2213 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 2214 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 2215 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2216 ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 2217 for (i=0; i<pn; i++) { 2218 row = i + rstart; 2219 nnz = bi[i+1] - bi[i]; 2220 Jptr = bj + bi[i]; 2221 ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,NULL,INSERT_VALUES);CHKERRQ(ierr); 2222 } 2223 ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2224 ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2225 merge->bi = bi; 2226 merge->bj = bj; 2227 merge->coi = coi; 2228 merge->coj = coj; 2229 merge->buf_ri = buf_ri; 2230 merge->buf_rj = buf_rj; 2231 merge->owners_co = owners_co; 2232 2233 /* attach the supporting struct to Cmpi for reuse */ 2234 c = (Mat_MPIAIJ*)Cmpi->data; 2235 2236 c->ap = ptap; 2237 ptap->api = NULL; 2238 ptap->apj = NULL; 2239 ptap->merge = merge; 2240 ptap->apa = NULL; 2241 ptap->destroy = Cmpi->ops->destroy; 2242 ptap->duplicate = Cmpi->ops->duplicate; 2243 2244 Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 2245 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 2246 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 2247 2248 *C = Cmpi; 2249 #if defined(PETSC_USE_INFO) 2250 if (bi[pn] != 0) { 2251 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 2252 ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 2253 } else { 2254 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 2255 } 2256 #endif 2257 PetscFunctionReturn(0); 2258 } 2259