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