1 2 /* 3 Defines projective product routines where A is a MPIAIJ matrix 4 C = P^T * A * P 5 */ 6 7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8 #include <../src/mat/utils/freespace.h> 9 #include <../src/mat/impls/aij/mpi/mpiaij.h> 10 #include <petscbt.h> 11 #include <petsctime.h> 12 13 /* #define PTAP_PROFILE */ 14 15 PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer) 16 { 17 PetscErrorCode ierr; 18 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 19 Mat_PtAPMPI *ptap=a->ptap; 20 PetscBool iascii; 21 PetscViewerFormat format; 22 23 PetscFunctionBegin; 24 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 25 if (iascii) { 26 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 27 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 28 if (ptap->algType == 0) { 29 ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr); 30 } else if (ptap->algType == 1) { 31 ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr); 32 } 33 } 34 } 35 ierr = (ptap->view)(A,viewer);CHKERRQ(ierr); 36 PetscFunctionReturn(0); 37 } 38 39 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 40 { 41 PetscErrorCode ierr; 42 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 43 Mat_PtAPMPI *ptap=a->ptap; 44 45 PetscFunctionBegin; 46 if (ptap) { 47 Mat_Merge_SeqsToMPI *merge=ptap->merge; 48 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 49 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 50 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 51 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 52 ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 53 ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 54 ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 55 if (ptap->AP_loc) { /* used by alg_rap */ 56 Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 57 ierr = PetscFree(ap->i);CHKERRQ(ierr); 58 ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 59 ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 60 } else { /* used by alg_ptap */ 61 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 62 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 63 } 64 ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 65 ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 66 if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 67 68 if (merge) { /* used by alg_ptap */ 69 ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 70 ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 71 ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 72 ierr = PetscFree(merge->bi);CHKERRQ(ierr); 73 ierr = PetscFree(merge->bj);CHKERRQ(ierr); 74 ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 75 ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 76 ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 77 ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 78 ierr = PetscFree(merge->coi);CHKERRQ(ierr); 79 ierr = PetscFree(merge->coj);CHKERRQ(ierr); 80 ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 81 ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 82 ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 83 } 84 85 ierr = ptap->destroy(A);CHKERRQ(ierr); 86 ierr = PetscFree(ptap);CHKERRQ(ierr); 87 } 88 PetscFunctionReturn(0); 89 } 90 91 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 92 { 93 PetscErrorCode ierr; 94 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 95 Mat_PtAPMPI *ptap = a->ptap; 96 97 PetscFunctionBegin; 98 ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 99 (*M)->ops->destroy = ptap->destroy; 100 (*M)->ops->duplicate = ptap->duplicate; 101 (*M)->ops->view = ptap->view; 102 PetscFunctionReturn(0); 103 } 104 105 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 106 { 107 PetscErrorCode ierr; 108 PetscBool flg; 109 MPI_Comm comm; 110 #if !defined(PETSC_HAVE_HYPRE) 111 const char *algTypes[2] = {"scalable","nonscalable"}; 112 PetscInt nalg=2; 113 #else 114 const char *algTypes[3] = {"scalable","nonscalable","hypre"}; 115 PetscInt nalg=3; 116 #endif 117 PetscInt pN=P->cmap->N,alg=1; /* set default algorithm */ 118 119 PetscFunctionBegin; 120 /* check if matrix local sizes are compatible */ 121 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 122 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,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); 123 if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(comm,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); 124 125 if (scall == MAT_INITIAL_MATRIX) { 126 /* pick an algorithm */ 127 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 128 PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 129 ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 130 ierr = PetscOptionsEnd();CHKERRQ(ierr); 131 132 if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */ 133 MatInfo Ainfo,Pinfo; 134 PetscInt nz_local; 135 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 136 137 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 138 ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 139 nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); 140 141 if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 142 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 143 144 if (alg_scalable) { 145 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 146 } 147 } 148 149 switch (alg) { 150 case 1: 151 /* do R=P^T locally, then C=R*A*P */ 152 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 153 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 154 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 155 break; 156 #if defined(PETSC_HAVE_HYPRE) 157 case 2: 158 /* Use boomerAMGBuildCoarseOperator */ 159 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 160 ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 161 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 162 break; 163 #endif 164 default: 165 /* do R=P^T locally, then C=R*A*P */ 166 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 167 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr); 168 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 169 break; 170 } 171 } 172 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 173 ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 174 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 175 PetscFunctionReturn(0); 176 } 177 178 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C) 179 { 180 PetscErrorCode ierr; 181 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 182 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 183 Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 184 Mat_PtAPMPI *ptap = c->ptap; 185 Mat AP_loc,C_loc,C_oth; 186 PetscInt i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz; 187 PetscScalar *apa; 188 const PetscInt *cols; 189 const PetscScalar *vals; 190 191 PetscFunctionBegin; 192 ierr = MatZeroEntries(C);CHKERRQ(ierr); 193 194 /* 1) get R = Pd^T,Ro = Po^T */ 195 if (ptap->reuse == MAT_REUSE_MATRIX) { 196 ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 197 ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 198 } 199 200 /* 2) get AP_loc */ 201 AP_loc = ptap->AP_loc; 202 ap = (Mat_SeqAIJ*)AP_loc->data; 203 204 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 205 /*-----------------------------------------------------*/ 206 if (ptap->reuse == MAT_REUSE_MATRIX) { 207 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 208 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 209 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 210 } 211 212 /* 2-2) compute numeric A_loc*P - dominating part */ 213 /* ---------------------------------------------- */ 214 /* get data from symbolic products */ 215 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 216 if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 217 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!"); 218 219 api = ap->i; 220 apj = ap->j; 221 for (i=0; i<am; i++) { 222 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 223 apnz = api[i+1] - api[i]; 224 apa = ap->a + api[i]; 225 ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr); 226 AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa); 227 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 228 } 229 230 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 231 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 232 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 233 234 C_loc = ptap->C_loc; 235 C_oth = ptap->C_oth; 236 237 /* add C_loc and Co to to C */ 238 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 239 240 /* C_loc -> C */ 241 cm = C_loc->rmap->N; 242 c_seq = (Mat_SeqAIJ*)C_loc->data; 243 cols = c_seq->j; 244 vals = c_seq->a; 245 246 /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 247 /* when there are no off-processor parts. */ 248 /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 249 /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 250 /* a table, and other, more complex stuff has to be done. */ 251 if (C->assembled) { 252 C->was_assembled = PETSC_TRUE; 253 C->assembled = PETSC_FALSE; 254 } 255 if (C->was_assembled) { 256 for (i=0; i<cm; i++) { 257 ncols = c_seq->i[i+1] - c_seq->i[i]; 258 row = rstart + i; 259 ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 260 cols += ncols; vals += ncols; 261 } 262 } else { 263 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 264 } 265 266 /* Co -> C, off-processor part */ 267 cm = C_oth->rmap->N; 268 c_seq = (Mat_SeqAIJ*)C_oth->data; 269 cols = c_seq->j; 270 vals = c_seq->a; 271 for (i=0; i<cm; i++) { 272 ncols = c_seq->i[i+1] - c_seq->i[i]; 273 row = p->garray[i]; 274 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 275 cols += ncols; vals += ncols; 276 } 277 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 278 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 279 280 ptap->reuse = MAT_REUSE_MATRIX; 281 PetscFunctionReturn(0); 282 } 283 284 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C) 285 { 286 PetscErrorCode ierr; 287 Mat_PtAPMPI *ptap; 288 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 289 MPI_Comm comm; 290 PetscMPIInt size,rank; 291 Mat Cmpi,P_loc,P_oth; 292 PetscFreeSpaceList free_space=NULL,current_space=NULL; 293 PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 294 PetscInt *lnk,i,k,pnz,row,nsend; 295 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 296 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 297 PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 298 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 299 MPI_Request *swaits,*rwaits; 300 MPI_Status *sstatus,rstatus; 301 PetscLayout rowmap; 302 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 303 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 304 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi; 305 Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 306 PetscScalar *apv; 307 PetscTable ta; 308 #if defined(PETSC_USE_INFO) 309 PetscReal apfill; 310 #endif 311 312 PetscFunctionBegin; 313 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 314 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 315 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 316 317 if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 318 319 /* create symbolic parallel matrix Cmpi */ 320 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 321 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 322 323 /* create struct Mat_PtAPMPI and attached it to C later */ 324 ierr = PetscNew(&ptap);CHKERRQ(ierr); 325 ptap->reuse = MAT_INITIAL_MATRIX; 326 ptap->algType = 0; 327 328 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 329 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr); 330 /* get P_loc by taking all local rows of P */ 331 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr); 332 333 ptap->P_loc = P_loc; 334 ptap->P_oth = P_oth; 335 336 /* (0) compute Rd = Pd^T, Ro = Po^T */ 337 /* --------------------------------- */ 338 ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 339 ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 340 341 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 342 /* ----------------------------------------------------------------- */ 343 p_loc = (Mat_SeqAIJ*)P_loc->data; 344 if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data; 345 346 /* create and initialize a linked list */ 347 ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 348 MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta); 349 MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta); 350 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 351 352 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 353 354 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 355 if (ao) { 356 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 357 } else { 358 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 359 } 360 current_space = free_space; 361 nspacedouble = 0; 362 363 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 364 api[0] = 0; 365 for (i=0; i<am; i++) { 366 /* diagonal portion: Ad[i,:]*P */ 367 ai = ad->i; pi = p_loc->i; 368 nzi = ai[i+1] - ai[i]; 369 aj = ad->j + ai[i]; 370 for (j=0; j<nzi; j++) { 371 row = aj[j]; 372 pnz = pi[row+1] - pi[row]; 373 Jptr = p_loc->j + pi[row]; 374 /* add non-zero cols of P into the sorted linked list lnk */ 375 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 376 } 377 /* off-diagonal portion: Ao[i,:]*P */ 378 if (ao) { 379 ai = ao->i; pi = p_oth->i; 380 nzi = ai[i+1] - ai[i]; 381 aj = ao->j + ai[i]; 382 for (j=0; j<nzi; j++) { 383 row = aj[j]; 384 pnz = pi[row+1] - pi[row]; 385 Jptr = p_oth->j + pi[row]; 386 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 387 } 388 } 389 apnz = lnk[0]; 390 api[i+1] = api[i] + apnz; 391 392 /* if free space is not available, double the total space in the list */ 393 if (current_space->local_remaining<apnz) { 394 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 395 nspacedouble++; 396 } 397 398 /* Copy data into free space, then initialize lnk */ 399 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 400 401 current_space->array += apnz; 402 current_space->local_used += apnz; 403 current_space->local_remaining -= apnz; 404 } 405 /* Allocate space for apj and apv, initialize apj, and */ 406 /* destroy list of free space and other temporary array(s) */ 407 ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 408 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 409 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 410 411 /* Create AP_loc for reuse */ 412 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 413 414 #if defined(PETSC_USE_INFO) 415 if (ao) { 416 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 417 } else { 418 apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 419 } 420 ptap->AP_loc->info.mallocs = nspacedouble; 421 ptap->AP_loc->info.fill_ratio_given = fill; 422 ptap->AP_loc->info.fill_ratio_needed = apfill; 423 424 if (api[am]) { 425 ierr = PetscInfo3(ptap->AP_loc,"Scalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 426 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 427 } else { 428 ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 429 } 430 #endif 431 432 /* (2-1) compute symbolic Co = Ro*AP_loc */ 433 /* ------------------------------------ */ 434 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 435 436 /* (3) send coj of C_oth to other processors */ 437 /* ------------------------------------------ */ 438 /* determine row ownership */ 439 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 440 rowmap->n = pn; 441 rowmap->bs = 1; 442 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 443 owners = rowmap->range; 444 445 /* determine the number of messages to send, their lengths */ 446 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 447 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 448 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 449 450 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 451 coi = c_oth->i; coj = c_oth->j; 452 con = ptap->C_oth->rmap->n; 453 proc = 0; 454 for (i=0; i<con; i++) { 455 while (prmap[i] >= owners[proc+1]) proc++; 456 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 457 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 458 } 459 460 len = 0; /* max length of buf_si[], see (4) */ 461 owners_co[0] = 0; 462 nsend = 0; 463 for (proc=0; proc<size; proc++) { 464 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 465 if (len_s[proc]) { 466 nsend++; 467 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 468 len += len_si[proc]; 469 } 470 } 471 472 /* determine the number and length of messages to receive for coi and coj */ 473 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 474 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 475 476 /* post the Irecv and Isend of coj */ 477 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 478 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 479 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 480 for (proc=0, k=0; proc<size; proc++) { 481 if (!len_s[proc]) continue; 482 i = owners_co[proc]; 483 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 484 k++; 485 } 486 487 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 488 /* ---------------------------------------- */ 489 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 490 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 491 492 /* receives coj are complete */ 493 for (i=0; i<nrecv; i++) { 494 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 495 } 496 ierr = PetscFree(rwaits);CHKERRQ(ierr); 497 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 498 499 /* add received column indices into ta to update Crmax */ 500 for (k=0; k<nrecv; k++) {/* k-th received message */ 501 Jptr = buf_rj[k]; 502 for (j=0; j<len_r[k]; j++) { 503 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 504 } 505 } 506 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 507 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 508 509 /* (4) send and recv coi */ 510 /*-----------------------*/ 511 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 512 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 513 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 514 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 515 for (proc=0,k=0; proc<size; proc++) { 516 if (!len_s[proc]) continue; 517 /* form outgoing message for i-structure: 518 buf_si[0]: nrows to be sent 519 [1:nrows]: row index (global) 520 [nrows+1:2*nrows+1]: i-structure index 521 */ 522 /*-------------------------------------------*/ 523 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 524 buf_si_i = buf_si + nrows+1; 525 buf_si[0] = nrows; 526 buf_si_i[0] = 0; 527 nrows = 0; 528 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 529 nzi = coi[i+1] - coi[i]; 530 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 531 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 532 nrows++; 533 } 534 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 535 k++; 536 buf_si += len_si[proc]; 537 } 538 for (i=0; i<nrecv; i++) { 539 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 540 } 541 ierr = PetscFree(rwaits);CHKERRQ(ierr); 542 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 543 544 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 545 ierr = PetscFree(len_ri);CHKERRQ(ierr); 546 ierr = PetscFree(swaits);CHKERRQ(ierr); 547 ierr = PetscFree(buf_s);CHKERRQ(ierr); 548 549 /* (5) compute the local portion of Cmpi */ 550 /* ------------------------------------------ */ 551 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 552 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 553 current_space = free_space; 554 555 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 556 for (k=0; k<nrecv; k++) { 557 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 558 nrows = *buf_ri_k[k]; 559 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 560 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 561 } 562 563 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 564 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 565 for (i=0; i<pn; i++) { 566 /* add C_loc into Cmpi */ 567 nzi = c_loc->i[i+1] - c_loc->i[i]; 568 Jptr = c_loc->j + c_loc->i[i]; 569 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 570 571 /* add received col data into lnk */ 572 for (k=0; k<nrecv; k++) { /* k-th received message */ 573 if (i == *nextrow[k]) { /* i-th row */ 574 nzi = *(nextci[k]+1) - *nextci[k]; 575 Jptr = buf_rj[k] + *nextci[k]; 576 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 577 nextrow[k]++; nextci[k]++; 578 } 579 } 580 nzi = lnk[0]; 581 582 /* copy data into free space, then initialize lnk */ 583 ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr); 584 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 585 } 586 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 587 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 588 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 589 590 /* local sizes and preallocation */ 591 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 592 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 593 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 594 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 595 596 /* members in merge */ 597 ierr = PetscFree(id_r);CHKERRQ(ierr); 598 ierr = PetscFree(len_r);CHKERRQ(ierr); 599 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 600 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 601 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 602 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 603 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 604 605 /* attach the supporting struct to Cmpi for reuse */ 606 c = (Mat_MPIAIJ*)Cmpi->data; 607 c->ptap = ptap; 608 ptap->duplicate = Cmpi->ops->duplicate; 609 ptap->destroy = Cmpi->ops->destroy; 610 ptap->view = Cmpi->ops->view; 611 612 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 613 Cmpi->assembled = PETSC_FALSE; 614 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 615 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 616 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 617 Cmpi->ops->view = MatView_MPIAIJ_PtAP; 618 *C = Cmpi; 619 PetscFunctionReturn(0); 620 } 621 622 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 623 { 624 PetscErrorCode ierr; 625 Mat_PtAPMPI *ptap; 626 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 627 MPI_Comm comm; 628 PetscMPIInt size,rank; 629 Mat Cmpi; 630 PetscFreeSpaceList free_space=NULL,current_space=NULL; 631 PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 632 PetscInt *lnk,i,k,pnz,row,nsend; 633 PetscBT lnkbt; 634 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 635 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 636 PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 637 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 638 MPI_Request *swaits,*rwaits; 639 MPI_Status *sstatus,rstatus; 640 PetscLayout rowmap; 641 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 642 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 643 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 644 Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 645 PetscScalar *apv; 646 PetscTable ta; 647 #if defined(PETSC_USE_INFO) 648 PetscReal apfill; 649 #endif 650 651 PetscFunctionBegin; 652 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 653 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 654 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 655 656 if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 657 658 /* create symbolic parallel matrix Cmpi */ 659 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 660 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 661 662 /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 663 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 664 665 /* create struct Mat_PtAPMPI and attached it to C later */ 666 ierr = PetscNew(&ptap);CHKERRQ(ierr); 667 ptap->reuse = MAT_INITIAL_MATRIX; 668 ptap->algType = 1; 669 670 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 671 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 672 /* get P_loc by taking all local rows of P */ 673 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 674 675 /* (0) compute Rd = Pd^T, Ro = Po^T */ 676 /* --------------------------------- */ 677 ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 678 ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 679 680 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 681 /* ----------------------------------------------------------------- */ 682 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 683 if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 684 685 /* create and initialize a linked list */ 686 ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 687 MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 688 MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 689 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 690 /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */ 691 692 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 693 694 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 695 if (ao) { 696 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 697 } else { 698 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 699 } 700 current_space = free_space; 701 nspacedouble = 0; 702 703 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 704 api[0] = 0; 705 for (i=0; i<am; i++) { 706 /* diagonal portion: Ad[i,:]*P */ 707 ai = ad->i; pi = p_loc->i; 708 nzi = ai[i+1] - ai[i]; 709 aj = ad->j + ai[i]; 710 for (j=0; j<nzi; j++) { 711 row = aj[j]; 712 pnz = pi[row+1] - pi[row]; 713 Jptr = p_loc->j + pi[row]; 714 /* add non-zero cols of P into the sorted linked list lnk */ 715 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 716 } 717 /* off-diagonal portion: Ao[i,:]*P */ 718 if (ao) { 719 ai = ao->i; pi = p_oth->i; 720 nzi = ai[i+1] - ai[i]; 721 aj = ao->j + ai[i]; 722 for (j=0; j<nzi; j++) { 723 row = aj[j]; 724 pnz = pi[row+1] - pi[row]; 725 Jptr = p_oth->j + pi[row]; 726 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 727 } 728 } 729 apnz = lnk[0]; 730 api[i+1] = api[i] + apnz; 731 if (ap_rmax < apnz) ap_rmax = apnz; 732 733 /* if free space is not available, double the total space in the list */ 734 if (current_space->local_remaining<apnz) { 735 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 736 nspacedouble++; 737 } 738 739 /* Copy data into free space, then initialize lnk */ 740 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 741 742 current_space->array += apnz; 743 current_space->local_used += apnz; 744 current_space->local_remaining -= apnz; 745 } 746 /* Allocate space for apj and apv, initialize apj, and */ 747 /* destroy list of free space and other temporary array(s) */ 748 ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 749 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 750 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 751 752 /* Create AP_loc for reuse */ 753 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 754 755 #if defined(PETSC_USE_INFO) 756 if (ao) { 757 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 758 } else { 759 apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 760 } 761 ptap->AP_loc->info.mallocs = nspacedouble; 762 ptap->AP_loc->info.fill_ratio_given = fill; 763 ptap->AP_loc->info.fill_ratio_needed = apfill; 764 765 if (api[am]) { 766 ierr = PetscInfo3(ptap->AP_loc,"Nonscalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 767 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 768 } else { 769 ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 770 } 771 #endif 772 773 /* (2-1) compute symbolic Co = Ro*AP_loc */ 774 /* ------------------------------------ */ 775 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 776 777 /* (3) send coj of C_oth to other processors */ 778 /* ------------------------------------------ */ 779 /* determine row ownership */ 780 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 781 rowmap->n = pn; 782 rowmap->bs = 1; 783 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 784 owners = rowmap->range; 785 786 /* determine the number of messages to send, their lengths */ 787 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 788 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 789 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 790 791 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 792 coi = c_oth->i; coj = c_oth->j; 793 con = ptap->C_oth->rmap->n; 794 proc = 0; 795 for (i=0; i<con; i++) { 796 while (prmap[i] >= owners[proc+1]) proc++; 797 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 798 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 799 } 800 801 len = 0; /* max length of buf_si[], see (4) */ 802 owners_co[0] = 0; 803 nsend = 0; 804 for (proc=0; proc<size; proc++) { 805 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 806 if (len_s[proc]) { 807 nsend++; 808 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 809 len += len_si[proc]; 810 } 811 } 812 813 /* determine the number and length of messages to receive for coi and coj */ 814 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 815 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 816 817 /* post the Irecv and Isend of coj */ 818 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 819 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 820 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 821 for (proc=0, k=0; proc<size; proc++) { 822 if (!len_s[proc]) continue; 823 i = owners_co[proc]; 824 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 825 k++; 826 } 827 828 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 829 /* ---------------------------------------- */ 830 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 831 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 832 833 /* receives coj are complete */ 834 for (i=0; i<nrecv; i++) { 835 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 836 } 837 ierr = PetscFree(rwaits);CHKERRQ(ierr); 838 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 839 840 /* add received column indices into ta to update Crmax */ 841 for (k=0; k<nrecv; k++) {/* k-th received message */ 842 Jptr = buf_rj[k]; 843 for (j=0; j<len_r[k]; j++) { 844 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 845 } 846 } 847 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 848 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 849 850 /* (4) send and recv coi */ 851 /*-----------------------*/ 852 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 853 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 854 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 855 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 856 for (proc=0,k=0; proc<size; proc++) { 857 if (!len_s[proc]) continue; 858 /* form outgoing message for i-structure: 859 buf_si[0]: nrows to be sent 860 [1:nrows]: row index (global) 861 [nrows+1:2*nrows+1]: i-structure index 862 */ 863 /*-------------------------------------------*/ 864 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 865 buf_si_i = buf_si + nrows+1; 866 buf_si[0] = nrows; 867 buf_si_i[0] = 0; 868 nrows = 0; 869 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 870 nzi = coi[i+1] - coi[i]; 871 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 872 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 873 nrows++; 874 } 875 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 876 k++; 877 buf_si += len_si[proc]; 878 } 879 for (i=0; i<nrecv; i++) { 880 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 881 } 882 ierr = PetscFree(rwaits);CHKERRQ(ierr); 883 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 884 885 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 886 ierr = PetscFree(len_ri);CHKERRQ(ierr); 887 ierr = PetscFree(swaits);CHKERRQ(ierr); 888 ierr = PetscFree(buf_s);CHKERRQ(ierr); 889 890 /* (5) compute the local portion of Cmpi */ 891 /* ------------------------------------------ */ 892 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 893 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 894 current_space = free_space; 895 896 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 897 for (k=0; k<nrecv; k++) { 898 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 899 nrows = *buf_ri_k[k]; 900 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 901 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 902 } 903 904 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 905 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 906 for (i=0; i<pn; i++) { 907 /* add C_loc into Cmpi */ 908 nzi = c_loc->i[i+1] - c_loc->i[i]; 909 Jptr = c_loc->j + c_loc->i[i]; 910 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 911 912 /* add received col data into lnk */ 913 for (k=0; k<nrecv; k++) { /* k-th received message */ 914 if (i == *nextrow[k]) { /* i-th row */ 915 nzi = *(nextci[k]+1) - *nextci[k]; 916 Jptr = buf_rj[k] + *nextci[k]; 917 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 918 nextrow[k]++; nextci[k]++; 919 } 920 } 921 nzi = lnk[0]; 922 923 /* copy data into free space, then initialize lnk */ 924 ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 925 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 926 } 927 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 928 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 929 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 930 931 /* local sizes and preallocation */ 932 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 933 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 934 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 935 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 936 937 /* members in merge */ 938 ierr = PetscFree(id_r);CHKERRQ(ierr); 939 ierr = PetscFree(len_r);CHKERRQ(ierr); 940 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 941 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 942 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 943 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 944 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 945 946 /* attach the supporting struct to Cmpi for reuse */ 947 c = (Mat_MPIAIJ*)Cmpi->data; 948 c->ptap = ptap; 949 ptap->duplicate = Cmpi->ops->duplicate; 950 ptap->destroy = Cmpi->ops->destroy; 951 ptap->view = Cmpi->ops->view; 952 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 953 954 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 955 Cmpi->assembled = PETSC_FALSE; 956 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 957 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 958 Cmpi->ops->view = MatView_MPIAIJ_PtAP; 959 *C = Cmpi; 960 PetscFunctionReturn(0); 961 } 962 963 964 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 965 { 966 PetscErrorCode ierr; 967 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 968 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 969 Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 970 Mat_PtAPMPI *ptap = c->ptap; 971 Mat AP_loc,C_loc,C_oth; 972 PetscInt i,rstart,rend,cm,ncols,row; 973 PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 974 PetscScalar *apa; 975 const PetscInt *cols; 976 const PetscScalar *vals; 977 978 PetscFunctionBegin; 979 ierr = MatZeroEntries(C);CHKERRQ(ierr); 980 /* 1) get R = Pd^T,Ro = Po^T */ 981 if (ptap->reuse == MAT_REUSE_MATRIX) { 982 ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 983 ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 984 } 985 986 /* 2) get AP_loc */ 987 AP_loc = ptap->AP_loc; 988 ap = (Mat_SeqAIJ*)AP_loc->data; 989 990 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 991 /*-----------------------------------------------------*/ 992 if (ptap->reuse == MAT_REUSE_MATRIX) { 993 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 994 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 995 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 996 } 997 998 /* 2-2) compute numeric A_loc*P - dominating part */ 999 /* ---------------------------------------------- */ 1000 /* get data from symbolic products */ 1001 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1002 if (ptap->P_oth) { 1003 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1004 } 1005 apa = ptap->apa; 1006 api = ap->i; 1007 apj = ap->j; 1008 for (i=0; i<am; i++) { 1009 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 1010 AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 1011 apnz = api[i+1] - api[i]; 1012 for (j=0; j<apnz; j++) { 1013 col = apj[j+api[i]]; 1014 ap->a[j+ap->i[i]] = apa[col]; 1015 apa[col] = 0.0; 1016 } 1017 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1018 } 1019 1020 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 1021 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 1022 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 1023 C_loc = ptap->C_loc; 1024 C_oth = ptap->C_oth; 1025 1026 /* add C_loc and Co to to C */ 1027 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 1028 1029 /* C_loc -> C */ 1030 cm = C_loc->rmap->N; 1031 c_seq = (Mat_SeqAIJ*)C_loc->data; 1032 cols = c_seq->j; 1033 vals = c_seq->a; 1034 1035 1036 /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 1037 /* when there are no off-processor parts. */ 1038 /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 1039 /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 1040 /* a table, and other, more complex stuff has to be done. */ 1041 if (C->assembled) { 1042 C->was_assembled = PETSC_TRUE; 1043 C->assembled = PETSC_FALSE; 1044 } 1045 if (C->was_assembled) { 1046 for (i=0; i<cm; i++) { 1047 ncols = c_seq->i[i+1] - c_seq->i[i]; 1048 row = rstart + i; 1049 ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1050 cols += ncols; vals += ncols; 1051 } 1052 } else { 1053 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 1054 } 1055 1056 /* Co -> C, off-processor part */ 1057 cm = C_oth->rmap->N; 1058 c_seq = (Mat_SeqAIJ*)C_oth->data; 1059 cols = c_seq->j; 1060 vals = c_seq->a; 1061 for (i=0; i<cm; i++) { 1062 ncols = c_seq->i[i+1] - c_seq->i[i]; 1063 row = p->garray[i]; 1064 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1065 cols += ncols; vals += ncols; 1066 } 1067 1068 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1069 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1070 1071 ptap->reuse = MAT_REUSE_MATRIX; 1072 PetscFunctionReturn(0); 1073 } 1074