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