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