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