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