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