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 MatDestroy_MPIAIJ_PtAP(Mat A) 16 { 17 PetscErrorCode ierr; 18 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 19 Mat_PtAPMPI *ptap=a->ptap; 20 21 PetscFunctionBegin; 22 if (ptap) { 23 Mat_Merge_SeqsToMPI *merge=ptap->merge; 24 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 25 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 26 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 27 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 28 ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 29 ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 30 ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 31 if (ptap->AP_loc) { /* used by alg_rap */ 32 Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 33 ierr = PetscFree(ap->i);CHKERRQ(ierr); 34 ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 35 ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 36 } else { /* used by alg_ptap */ 37 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 38 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 39 } 40 ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 41 ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 42 if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 43 44 if (merge) { /* used by alg_ptap */ 45 ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 46 ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 47 ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 48 ierr = PetscFree(merge->bi);CHKERRQ(ierr); 49 ierr = PetscFree(merge->bj);CHKERRQ(ierr); 50 ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 51 ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 52 ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 53 ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 54 ierr = PetscFree(merge->coi);CHKERRQ(ierr); 55 ierr = PetscFree(merge->coj);CHKERRQ(ierr); 56 ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 57 ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 58 ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 59 } 60 61 ierr = ptap->destroy(A);CHKERRQ(ierr); 62 ierr = PetscFree(ptap);CHKERRQ(ierr); 63 } 64 PetscFunctionReturn(0); 65 } 66 67 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 68 { 69 PetscErrorCode ierr; 70 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 71 Mat_PtAPMPI *ptap = a->ptap; 72 73 PetscFunctionBegin; 74 ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 75 (*M)->ops->destroy = ptap->destroy; 76 (*M)->ops->duplicate = ptap->duplicate; 77 PetscFunctionReturn(0); 78 } 79 80 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 81 { 82 PetscErrorCode ierr; 83 PetscBool rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */ 84 MPI_Comm comm; 85 86 PetscFunctionBegin; 87 /* check if matrix local sizes are compatible */ 88 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 89 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); 90 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); 91 92 ierr = PetscOptionsGetBool(NULL,NULL,"-matrap",&rap,NULL);CHKERRQ(ierr); 93 if (scall == MAT_INITIAL_MATRIX) { 94 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 95 if (rap) { /* do R=P^T locally, then C=R*A*P */ 96 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 97 } else { /* do P^T*A*P */ 98 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr); 99 } 100 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 101 } 102 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 103 ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 104 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 105 PetscFunctionReturn(0); 106 } 107 108 #if defined(PETSC_HAVE_HYPRE) 109 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 110 #endif 111 112 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C) 113 { 114 PetscErrorCode ierr; 115 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 116 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 117 Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 118 Mat_PtAPMPI *ptap = c->ptap; 119 Mat AP_loc,C_loc,C_oth; 120 PetscInt i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz; 121 PetscScalar *apa; 122 const PetscInt *cols; 123 const PetscScalar *vals; 124 125 PetscFunctionBegin; 126 ierr = MatZeroEntries(C);CHKERRQ(ierr); 127 128 /* 1) get R = Pd^T,Ro = Po^T */ 129 if (ptap->reuse == MAT_REUSE_MATRIX) { 130 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 131 ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 132 } 133 134 /* 2) get AP_loc */ 135 AP_loc = ptap->AP_loc; 136 ap = (Mat_SeqAIJ*)AP_loc->data; 137 138 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 139 /*-----------------------------------------------------*/ 140 if (ptap->reuse == MAT_REUSE_MATRIX) { 141 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 142 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 143 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 144 } 145 146 /* 2-2) compute numeric A_loc*P - dominating part */ 147 /* ---------------------------------------------- */ 148 /* get data from symbolic products */ 149 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 150 if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 151 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!"); 152 153 api = ap->i; 154 apj = ap->j; 155 for (i=0; i<am; i++) { 156 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 157 apnz = api[i+1] - api[i]; 158 apa = ap->a + api[i]; 159 ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr); 160 AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa); 161 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 162 } 163 164 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 165 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 166 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 167 168 C_loc = ptap->C_loc; 169 C_oth = ptap->C_oth; 170 171 /* add C_loc and Co to to C */ 172 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 173 174 /* C_loc -> C */ 175 cm = C_loc->rmap->N; 176 c_seq = (Mat_SeqAIJ*)C_loc->data; 177 cols = c_seq->j; 178 vals = c_seq->a; 179 for (i=0; i<cm; i++) { 180 ncols = c_seq->i[i+1] - c_seq->i[i]; 181 row = rstart + i; 182 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 183 cols += ncols; vals += ncols; 184 } 185 186 /* Co -> C, off-processor part */ 187 cm = C_oth->rmap->N; 188 c_seq = (Mat_SeqAIJ*)C_oth->data; 189 cols = c_seq->j; 190 vals = c_seq->a; 191 for (i=0; i<cm; i++) { 192 ncols = c_seq->i[i+1] - c_seq->i[i]; 193 row = p->garray[i]; 194 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 195 cols += ncols; vals += ncols; 196 } 197 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 198 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 199 200 ptap->reuse = MAT_REUSE_MATRIX; 201 PetscFunctionReturn(0); 202 } 203 204 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C) 205 { 206 PetscErrorCode ierr; 207 Mat_PtAPMPI *ptap; 208 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 209 MPI_Comm comm; 210 PetscMPIInt size,rank; 211 Mat Cmpi,P_loc,P_oth; 212 PetscFreeSpaceList free_space=NULL,current_space=NULL; 213 PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 214 PetscInt *lnk,i,k,pnz,row,nsend; 215 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 216 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 217 PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 218 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 219 MPI_Request *swaits,*rwaits; 220 MPI_Status *sstatus,rstatus; 221 PetscLayout rowmap; 222 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 223 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 224 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi; 225 Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 226 PetscScalar *apv; 227 PetscTable ta; 228 #if defined(PETSC_USE_INFO) 229 PetscReal apfill; 230 #endif 231 232 PetscFunctionBegin; 233 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 234 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 235 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 236 237 if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 238 239 /* create symbolic parallel matrix Cmpi */ 240 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 241 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 242 243 /* create struct Mat_PtAPMPI and attached it to C later */ 244 ierr = PetscNew(&ptap);CHKERRQ(ierr); 245 ptap->reuse = MAT_INITIAL_MATRIX; 246 247 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 248 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr); 249 /* get P_loc by taking all local rows of P */ 250 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr); 251 252 ptap->P_loc = P_loc; 253 ptap->P_oth = P_oth; 254 255 /* (0) compute Rd = Pd^T, Ro = Po^T */ 256 /* --------------------------------- */ 257 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 258 ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 259 260 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 261 /* ----------------------------------------------------------------- */ 262 p_loc = (Mat_SeqAIJ*)P_loc->data; 263 if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data; 264 265 /* create and initialize a linked list */ 266 ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 267 MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta); 268 MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta); 269 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 270 271 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 272 273 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 274 if (ao) { 275 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 276 } else { 277 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 278 } 279 current_space = free_space; 280 nspacedouble = 0; 281 282 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 283 api[0] = 0; 284 for (i=0; i<am; i++) { 285 /* diagonal portion: Ad[i,:]*P */ 286 ai = ad->i; pi = p_loc->i; 287 nzi = ai[i+1] - ai[i]; 288 aj = ad->j + ai[i]; 289 for (j=0; j<nzi; j++) { 290 row = aj[j]; 291 pnz = pi[row+1] - pi[row]; 292 Jptr = p_loc->j + pi[row]; 293 /* add non-zero cols of P into the sorted linked list lnk */ 294 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 295 } 296 /* off-diagonal portion: Ao[i,:]*P */ 297 if (ao) { 298 ai = ao->i; pi = p_oth->i; 299 nzi = ai[i+1] - ai[i]; 300 aj = ao->j + ai[i]; 301 for (j=0; j<nzi; j++) { 302 row = aj[j]; 303 pnz = pi[row+1] - pi[row]; 304 Jptr = p_oth->j + pi[row]; 305 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 306 } 307 } 308 apnz = lnk[0]; 309 api[i+1] = api[i] + apnz; 310 311 /* if free space is not available, double the total space in the list */ 312 if (current_space->local_remaining<apnz) { 313 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 314 nspacedouble++; 315 } 316 317 /* Copy data into free space, then initialize lnk */ 318 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 319 320 current_space->array += apnz; 321 current_space->local_used += apnz; 322 current_space->local_remaining -= apnz; 323 } 324 /* Allocate space for apj and apv, initialize apj, and */ 325 /* destroy list of free space and other temporary array(s) */ 326 ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 327 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 328 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 329 330 /* Create AP_loc for reuse */ 331 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 332 333 #if defined(PETSC_USE_INFO) 334 if (ao) { 335 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 336 } else { 337 apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 338 } 339 ptap->AP_loc->info.mallocs = nspacedouble; 340 ptap->AP_loc->info.fill_ratio_given = fill; 341 ptap->AP_loc->info.fill_ratio_needed = apfill; 342 343 if (api[am]) { 344 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); 345 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 346 } else { 347 ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 348 } 349 #endif 350 351 /* (2-1) compute symbolic Co = Ro*AP_loc */ 352 /* ------------------------------------ */ 353 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 354 355 /* (3) send coj of C_oth to other processors */ 356 /* ------------------------------------------ */ 357 /* determine row ownership */ 358 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 359 rowmap->n = pn; 360 rowmap->bs = 1; 361 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 362 owners = rowmap->range; 363 364 /* determine the number of messages to send, their lengths */ 365 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 366 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 367 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 368 369 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 370 coi = c_oth->i; coj = c_oth->j; 371 con = ptap->C_oth->rmap->n; 372 proc = 0; 373 for (i=0; i<con; i++) { 374 while (prmap[i] >= owners[proc+1]) proc++; 375 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 376 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 377 } 378 379 len = 0; /* max length of buf_si[], see (4) */ 380 owners_co[0] = 0; 381 nsend = 0; 382 for (proc=0; proc<size; proc++) { 383 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 384 if (len_s[proc]) { 385 nsend++; 386 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 387 len += len_si[proc]; 388 } 389 } 390 391 /* determine the number and length of messages to receive for coi and coj */ 392 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 393 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 394 395 /* post the Irecv and Isend of coj */ 396 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 397 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 398 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 399 for (proc=0, k=0; proc<size; proc++) { 400 if (!len_s[proc]) continue; 401 i = owners_co[proc]; 402 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 403 k++; 404 } 405 406 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 407 /* ---------------------------------------- */ 408 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 409 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 410 411 /* receives coj are complete */ 412 for (i=0; i<nrecv; i++) { 413 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 414 } 415 ierr = PetscFree(rwaits);CHKERRQ(ierr); 416 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 417 418 /* add received column indices into ta to update Crmax */ 419 for (k=0; k<nrecv; k++) {/* k-th received message */ 420 Jptr = buf_rj[k]; 421 for (j=0; j<len_r[k]; j++) { 422 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 423 } 424 } 425 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 426 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 427 428 /* (4) send and recv coi */ 429 /*-----------------------*/ 430 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 431 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 432 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 433 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 434 for (proc=0,k=0; proc<size; proc++) { 435 if (!len_s[proc]) continue; 436 /* form outgoing message for i-structure: 437 buf_si[0]: nrows to be sent 438 [1:nrows]: row index (global) 439 [nrows+1:2*nrows+1]: i-structure index 440 */ 441 /*-------------------------------------------*/ 442 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 443 buf_si_i = buf_si + nrows+1; 444 buf_si[0] = nrows; 445 buf_si_i[0] = 0; 446 nrows = 0; 447 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 448 nzi = coi[i+1] - coi[i]; 449 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 450 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 451 nrows++; 452 } 453 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 454 k++; 455 buf_si += len_si[proc]; 456 } 457 for (i=0; i<nrecv; i++) { 458 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 459 } 460 ierr = PetscFree(rwaits);CHKERRQ(ierr); 461 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 462 463 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 464 ierr = PetscFree(len_ri);CHKERRQ(ierr); 465 ierr = PetscFree(swaits);CHKERRQ(ierr); 466 ierr = PetscFree(buf_s);CHKERRQ(ierr); 467 468 /* (5) compute the local portion of Cmpi */ 469 /* ------------------------------------------ */ 470 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 471 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 472 current_space = free_space; 473 474 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 475 for (k=0; k<nrecv; k++) { 476 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 477 nrows = *buf_ri_k[k]; 478 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 479 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 480 } 481 482 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 483 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 484 for (i=0; i<pn; i++) { 485 /* add C_loc into Cmpi */ 486 nzi = c_loc->i[i+1] - c_loc->i[i]; 487 Jptr = c_loc->j + c_loc->i[i]; 488 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 489 490 /* add received col data into lnk */ 491 for (k=0; k<nrecv; k++) { /* k-th received message */ 492 if (i == *nextrow[k]) { /* i-th row */ 493 nzi = *(nextci[k]+1) - *nextci[k]; 494 Jptr = buf_rj[k] + *nextci[k]; 495 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 496 nextrow[k]++; nextci[k]++; 497 } 498 } 499 nzi = lnk[0]; 500 501 /* copy data into free space, then initialize lnk */ 502 ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr); 503 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 504 } 505 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 506 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 507 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 508 509 /* local sizes and preallocation */ 510 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 511 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 512 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 513 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 514 515 /* members in merge */ 516 ierr = PetscFree(id_r);CHKERRQ(ierr); 517 ierr = PetscFree(len_r);CHKERRQ(ierr); 518 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 519 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 520 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 521 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 522 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 523 524 /* attach the supporting struct to Cmpi for reuse */ 525 c = (Mat_MPIAIJ*)Cmpi->data; 526 c->ptap = ptap; 527 ptap->duplicate = Cmpi->ops->duplicate; 528 ptap->destroy = Cmpi->ops->destroy; 529 530 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 531 Cmpi->assembled = PETSC_FALSE; 532 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 533 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 534 *C = Cmpi; 535 PetscFunctionReturn(0); 536 } 537 538 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 539 { 540 PetscErrorCode ierr; 541 Mat_PtAPMPI *ptap; 542 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 543 MPI_Comm comm; 544 PetscMPIInt size,rank; 545 Mat Cmpi; 546 PetscFreeSpaceList free_space=NULL,current_space=NULL; 547 PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 548 PetscInt *lnk,i,k,pnz,row,nsend; 549 PetscBT lnkbt; 550 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 551 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 552 PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 553 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 554 MPI_Request *swaits,*rwaits; 555 MPI_Status *sstatus,rstatus; 556 PetscLayout rowmap; 557 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 558 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 559 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 560 Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 561 PetscScalar *apv; 562 PetscTable ta; 563 #if defined(PETSC_HAVE_HYPRE) 564 const char *algTypes[3] = {"scalable","nonscalable","hypre"}; 565 PetscInt nalg = 3; 566 #else 567 const char *algTypes[2] = {"scalable","nonscalable"}; 568 PetscInt nalg = 2; 569 #endif 570 PetscInt alg = 1; /* set default algorithm */ 571 #if defined(PETSC_USE_INFO) 572 PetscReal apfill; 573 #endif 574 #if defined(PTAP_PROFILE) 575 PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 576 #endif 577 PetscBool flg; 578 579 PetscFunctionBegin; 580 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 581 582 /* pick an algorithm */ 583 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 584 PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 585 ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr); 586 ierr = PetscOptionsEnd();CHKERRQ(ierr); 587 588 if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */ 589 MatInfo Ainfo,Pinfo; 590 PetscInt nz_local; 591 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 592 593 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 594 ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 595 nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); 596 597 if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 598 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 599 600 if (alg_scalable) { 601 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 602 } 603 } 604 605 if (alg == 0) { 606 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr); 607 (*C)->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 608 PetscFunctionReturn(0); 609 610 #if defined(PETSC_HAVE_HYPRE) 611 } else if (alg == 2) { 612 /* Use boomerAMGBuildCoarseOperator */ 613 ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 614 PetscFunctionReturn(0); 615 #endif 616 } 617 618 #if defined(PTAP_PROFILE) 619 ierr = PetscTime(&t0);CHKERRQ(ierr); 620 #endif 621 622 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 623 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 624 625 if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 626 627 /* create symbolic parallel matrix Cmpi */ 628 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 629 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 630 631 /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 632 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 633 634 /* create struct Mat_PtAPMPI and attached it to C later */ 635 ierr = PetscNew(&ptap);CHKERRQ(ierr); 636 ptap->reuse = MAT_INITIAL_MATRIX; 637 638 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 639 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 640 /* get P_loc by taking all local rows of P */ 641 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 642 643 /* (0) compute Rd = Pd^T, Ro = Po^T */ 644 /* --------------------------------- */ 645 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 646 ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 647 #if defined(PTAP_PROFILE) 648 ierr = PetscTime(&t11);CHKERRQ(ierr); 649 #endif 650 651 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 652 /* ----------------------------------------------------------------- */ 653 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 654 if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 655 656 /* create and initialize a linked list */ 657 ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 658 MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 659 MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 660 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 661 /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */ 662 663 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 664 665 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 666 if (ao) { 667 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 668 } else { 669 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 670 } 671 current_space = free_space; 672 nspacedouble = 0; 673 674 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 675 api[0] = 0; 676 for (i=0; i<am; i++) { 677 /* diagonal portion: Ad[i,:]*P */ 678 ai = ad->i; pi = p_loc->i; 679 nzi = ai[i+1] - ai[i]; 680 aj = ad->j + ai[i]; 681 for (j=0; j<nzi; j++) { 682 row = aj[j]; 683 pnz = pi[row+1] - pi[row]; 684 Jptr = p_loc->j + pi[row]; 685 /* add non-zero cols of P into the sorted linked list lnk */ 686 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 687 } 688 /* off-diagonal portion: Ao[i,:]*P */ 689 if (ao) { 690 ai = ao->i; pi = p_oth->i; 691 nzi = ai[i+1] - ai[i]; 692 aj = ao->j + ai[i]; 693 for (j=0; j<nzi; j++) { 694 row = aj[j]; 695 pnz = pi[row+1] - pi[row]; 696 Jptr = p_oth->j + pi[row]; 697 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 698 } 699 } 700 apnz = lnk[0]; 701 api[i+1] = api[i] + apnz; 702 if (ap_rmax < apnz) ap_rmax = apnz; 703 704 /* if free space is not available, double the total space in the list */ 705 if (current_space->local_remaining<apnz) { 706 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 707 nspacedouble++; 708 } 709 710 /* Copy data into free space, then initialize lnk */ 711 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 712 713 current_space->array += apnz; 714 current_space->local_used += apnz; 715 current_space->local_remaining -= apnz; 716 } 717 /* Allocate space for apj and apv, initialize apj, and */ 718 /* destroy list of free space and other temporary array(s) */ 719 ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 720 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 721 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 722 723 /* Create AP_loc for reuse */ 724 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 725 726 #if defined(PETSC_USE_INFO) 727 if (ao) { 728 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 729 } else { 730 apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 731 } 732 ptap->AP_loc->info.mallocs = nspacedouble; 733 ptap->AP_loc->info.fill_ratio_given = fill; 734 ptap->AP_loc->info.fill_ratio_needed = apfill; 735 736 if (api[am]) { 737 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); 738 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 739 } else { 740 ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 741 } 742 #endif 743 744 #if defined(PTAP_PROFILE) 745 ierr = PetscTime(&t12);CHKERRQ(ierr); 746 #endif 747 748 /* (2-1) compute symbolic Co = Ro*AP_loc */ 749 /* ------------------------------------ */ 750 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 751 #if defined(PTAP_PROFILE) 752 ierr = PetscTime(&t1);CHKERRQ(ierr); 753 #endif 754 755 /* (3) send coj of C_oth to other processors */ 756 /* ------------------------------------------ */ 757 /* determine row ownership */ 758 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 759 rowmap->n = pn; 760 rowmap->bs = 1; 761 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 762 owners = rowmap->range; 763 764 /* determine the number of messages to send, their lengths */ 765 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 766 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 767 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 768 769 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 770 coi = c_oth->i; coj = c_oth->j; 771 con = ptap->C_oth->rmap->n; 772 proc = 0; 773 for (i=0; i<con; i++) { 774 while (prmap[i] >= owners[proc+1]) proc++; 775 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 776 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 777 } 778 779 len = 0; /* max length of buf_si[], see (4) */ 780 owners_co[0] = 0; 781 nsend = 0; 782 for (proc=0; proc<size; proc++) { 783 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 784 if (len_s[proc]) { 785 nsend++; 786 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 787 len += len_si[proc]; 788 } 789 } 790 791 /* determine the number and length of messages to receive for coi and coj */ 792 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 793 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 794 795 /* post the Irecv and Isend of coj */ 796 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 797 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 798 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 799 for (proc=0, k=0; proc<size; proc++) { 800 if (!len_s[proc]) continue; 801 i = owners_co[proc]; 802 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 803 k++; 804 } 805 806 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 807 /* ---------------------------------------- */ 808 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 809 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 810 811 /* receives coj are complete */ 812 for (i=0; i<nrecv; i++) { 813 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 814 } 815 ierr = PetscFree(rwaits);CHKERRQ(ierr); 816 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 817 818 /* add received column indices into ta to update Crmax */ 819 for (k=0; k<nrecv; k++) {/* k-th received message */ 820 Jptr = buf_rj[k]; 821 for (j=0; j<len_r[k]; j++) { 822 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 823 } 824 } 825 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 826 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 827 828 /* (4) send and recv coi */ 829 /*-----------------------*/ 830 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 831 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 832 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 833 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 834 for (proc=0,k=0; proc<size; proc++) { 835 if (!len_s[proc]) continue; 836 /* form outgoing message for i-structure: 837 buf_si[0]: nrows to be sent 838 [1:nrows]: row index (global) 839 [nrows+1:2*nrows+1]: i-structure index 840 */ 841 /*-------------------------------------------*/ 842 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 843 buf_si_i = buf_si + nrows+1; 844 buf_si[0] = nrows; 845 buf_si_i[0] = 0; 846 nrows = 0; 847 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 848 nzi = coi[i+1] - coi[i]; 849 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 850 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 851 nrows++; 852 } 853 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 854 k++; 855 buf_si += len_si[proc]; 856 } 857 for (i=0; i<nrecv; i++) { 858 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 859 } 860 ierr = PetscFree(rwaits);CHKERRQ(ierr); 861 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 862 863 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 864 ierr = PetscFree(len_ri);CHKERRQ(ierr); 865 ierr = PetscFree(swaits);CHKERRQ(ierr); 866 ierr = PetscFree(buf_s);CHKERRQ(ierr); 867 #if defined(PTAP_PROFILE) 868 ierr = PetscTime(&t2);CHKERRQ(ierr); 869 #endif 870 /* (5) compute the local portion of Cmpi */ 871 /* ------------------------------------------ */ 872 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 873 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 874 current_space = free_space; 875 876 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 877 for (k=0; k<nrecv; k++) { 878 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 879 nrows = *buf_ri_k[k]; 880 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 881 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 882 } 883 884 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 885 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 886 for (i=0; i<pn; i++) { 887 /* add C_loc into Cmpi */ 888 nzi = c_loc->i[i+1] - c_loc->i[i]; 889 Jptr = c_loc->j + c_loc->i[i]; 890 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 891 892 /* add received col data into lnk */ 893 for (k=0; k<nrecv; k++) { /* k-th received message */ 894 if (i == *nextrow[k]) { /* i-th row */ 895 nzi = *(nextci[k]+1) - *nextci[k]; 896 Jptr = buf_rj[k] + *nextci[k]; 897 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 898 nextrow[k]++; nextci[k]++; 899 } 900 } 901 nzi = lnk[0]; 902 903 /* copy data into free space, then initialize lnk */ 904 ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 905 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 906 } 907 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 908 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 909 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 910 #if defined(PTAP_PROFILE) 911 ierr = PetscTime(&t3);CHKERRQ(ierr); 912 #endif 913 914 /* local sizes and preallocation */ 915 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 916 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 917 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 918 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 919 920 /* members in merge */ 921 ierr = PetscFree(id_r);CHKERRQ(ierr); 922 ierr = PetscFree(len_r);CHKERRQ(ierr); 923 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 924 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 925 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 926 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 927 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 928 929 /* attach the supporting struct to Cmpi for reuse */ 930 c = (Mat_MPIAIJ*)Cmpi->data; 931 c->ptap = ptap; 932 ptap->duplicate = Cmpi->ops->duplicate; 933 ptap->destroy = Cmpi->ops->destroy; 934 935 if (alg == 1) { 936 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 937 } 938 939 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 940 Cmpi->assembled = PETSC_FALSE; 941 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 942 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 943 *C = Cmpi; 944 945 #if defined(PTAP_PROFILE) 946 ierr = PetscTime(&t4);CHKERRQ(ierr); 947 if (rank == 1) { 948 printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0); 949 } 950 #endif 951 PetscFunctionReturn(0); 952 } 953 954 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 955 { 956 PetscErrorCode ierr; 957 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 958 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 959 Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 960 Mat_PtAPMPI *ptap = c->ptap; 961 Mat AP_loc,C_loc,C_oth; 962 PetscInt i,rstart,rend,cm,ncols,row; 963 PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 964 PetscScalar *apa; 965 const PetscInt *cols; 966 const PetscScalar *vals; 967 #if defined(PTAP_PROFILE) 968 PetscMPIInt rank; 969 MPI_Comm comm; 970 PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 971 #endif 972 973 PetscFunctionBegin; 974 ierr = MatZeroEntries(C);CHKERRQ(ierr); 975 976 /* 1) get R = Pd^T,Ro = Po^T */ 977 #if defined(PTAP_PROFILE) 978 ierr = PetscTime(&t0);CHKERRQ(ierr); 979 #endif 980 if (ptap->reuse == MAT_REUSE_MATRIX) { 981 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 982 ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 983 } 984 #if defined(PTAP_PROFILE) 985 ierr = PetscTime(&t1);CHKERRQ(ierr); 986 eR = t1 - t0; 987 #endif 988 989 /* 2) get AP_loc */ 990 AP_loc = ptap->AP_loc; 991 ap = (Mat_SeqAIJ*)AP_loc->data; 992 993 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 994 /*-----------------------------------------------------*/ 995 if (ptap->reuse == MAT_REUSE_MATRIX) { 996 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 997 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 998 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 999 } 1000 1001 /* 2-2) compute numeric A_loc*P - dominating part */ 1002 /* ---------------------------------------------- */ 1003 /* get data from symbolic products */ 1004 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1005 if (ptap->P_oth) { 1006 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1007 } 1008 apa = ptap->apa; 1009 api = ap->i; 1010 apj = ap->j; 1011 for (i=0; i<am; i++) { 1012 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 1013 AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 1014 apnz = api[i+1] - api[i]; 1015 for (j=0; j<apnz; j++) { 1016 col = apj[j+api[i]]; 1017 ap->a[j+ap->i[i]] = apa[col]; 1018 apa[col] = 0.0; 1019 } 1020 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1021 } 1022 #if defined(PTAP_PROFILE) 1023 ierr = PetscTime(&t2);CHKERRQ(ierr); 1024 eAP = t2 - t1; 1025 #endif 1026 1027 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 1028 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 1029 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 1030 C_loc = ptap->C_loc; 1031 C_oth = ptap->C_oth; 1032 #if defined(PTAP_PROFILE) 1033 ierr = PetscTime(&t3);CHKERRQ(ierr); 1034 eCseq = t3 - t2; 1035 #endif 1036 1037 /* add C_loc and Co to to C */ 1038 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 1039 1040 /* C_loc -> C */ 1041 cm = C_loc->rmap->N; 1042 c_seq = (Mat_SeqAIJ*)C_loc->data; 1043 cols = c_seq->j; 1044 vals = c_seq->a; 1045 for (i=0; i<cm; i++) { 1046 ncols = c_seq->i[i+1] - c_seq->i[i]; 1047 row = rstart + i; 1048 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1049 cols += ncols; vals += ncols; 1050 } 1051 1052 /* Co -> C, off-processor part */ 1053 cm = C_oth->rmap->N; 1054 c_seq = (Mat_SeqAIJ*)C_oth->data; 1055 cols = c_seq->j; 1056 vals = c_seq->a; 1057 for (i=0; i<cm; i++) { 1058 ncols = c_seq->i[i+1] - c_seq->i[i]; 1059 row = p->garray[i]; 1060 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1061 cols += ncols; vals += ncols; 1062 } 1063 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1064 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1065 #if defined(PTAP_PROFILE) 1066 ierr = PetscTime(&t4);CHKERRQ(ierr); 1067 eCmpi = t4 - t3; 1068 1069 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1070 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1071 if (rank==1) { 1072 ierr = PetscPrintf(MPI_COMM_SELF," R %g, AP %g, Cseq %g, Cmpi %g = %g\n", eR,eAP,eCseq,eCmpi,eR+eAP+eCseq+eCmpi);CHKERRQ(ierr); 1073 } 1074 #endif 1075 ptap->reuse = MAT_REUSE_MATRIX; 1076 PetscFunctionReturn(0); 1077 } 1078 1079 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C) 1080 { 1081 PetscErrorCode ierr; 1082 Mat Cmpi; 1083 Mat_PtAPMPI *ptap; 1084 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1085 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 1086 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1087 Mat_SeqAIJ *p_loc,*p_oth; 1088 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 1089 PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 1090 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1091 PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 1092 PetscBT lnkbt; 1093 MPI_Comm comm; 1094 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 1095 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1096 PetscInt len,proc,*dnz,*onz,*owners; 1097 PetscInt nzi,*pti,*ptj; 1098 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1099 MPI_Request *swaits,*rwaits; 1100 MPI_Status *sstatus,rstatus; 1101 Mat_Merge_SeqsToMPI *merge; 1102 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 1103 PetscReal afill=1.0,afill_tmp; 1104 1105 PetscFunctionBegin; 1106 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1107 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1108 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1109 1110 /* create struct Mat_PtAPMPI and attached it to C later */ 1111 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1112 ierr = PetscNew(&merge);CHKERRQ(ierr); 1113 ptap->merge = merge; 1114 ptap->reuse = MAT_INITIAL_MATRIX; 1115 1116 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1117 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1118 1119 /* get P_loc by taking all local rows of P */ 1120 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1121 1122 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1123 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1124 pi_loc = p_loc->i; pj_loc = p_loc->j; 1125 pi_oth = p_oth->i; pj_oth = p_oth->j; 1126 1127 /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 1128 /*--------------------------------------------------------------------------*/ 1129 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 1130 api[0] = 0; 1131 1132 /* create and initialize a linked list */ 1133 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1134 1135 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 1136 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 1137 current_space = free_space; 1138 1139 for (i=0; i<am; i++) { 1140 /* diagonal portion of A */ 1141 nzi = adi[i+1] - adi[i]; 1142 aj = ad->j + adi[i]; 1143 for (j=0; j<nzi; j++) { 1144 row = aj[j]; 1145 pnz = pi_loc[row+1] - pi_loc[row]; 1146 Jptr = pj_loc + pi_loc[row]; 1147 /* add non-zero cols of P into the sorted linked list lnk */ 1148 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1149 } 1150 /* off-diagonal portion of A */ 1151 nzi = aoi[i+1] - aoi[i]; 1152 aj = ao->j + aoi[i]; 1153 for (j=0; j<nzi; j++) { 1154 row = aj[j]; 1155 pnz = pi_oth[row+1] - pi_oth[row]; 1156 Jptr = pj_oth + pi_oth[row]; 1157 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1158 } 1159 apnz = lnk[0]; 1160 api[i+1] = api[i] + apnz; 1161 if (ap_rmax < apnz) ap_rmax = apnz; 1162 1163 /* if free space is not available, double the total space in the list */ 1164 if (current_space->local_remaining<apnz) { 1165 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1166 nspacedouble++; 1167 } 1168 1169 /* Copy data into free space, then initialize lnk */ 1170 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1171 1172 current_space->array += apnz; 1173 current_space->local_used += apnz; 1174 current_space->local_remaining -= apnz; 1175 } 1176 1177 /* Allocate space for apj, initialize apj, and */ 1178 /* destroy list of free space and other temporary array(s) */ 1179 ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 1180 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 1181 afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 1182 if (afill_tmp > afill) afill = afill_tmp; 1183 1184 /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 1185 /*-----------------------------------------------------------------*/ 1186 ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 1187 1188 /* then, compute symbolic Co = (p->B)^T*AP */ 1189 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1190 >= (num of nonzero rows of C_seq) - pn */ 1191 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1192 coi[0] = 0; 1193 1194 /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 1195 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am])); 1196 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1197 current_space = free_space; 1198 1199 for (i=0; i<pon; i++) { 1200 pnz = poti[i+1] - poti[i]; 1201 ptJ = potj + poti[i]; 1202 for (j=0; j<pnz; j++) { 1203 row = ptJ[j]; /* row of AP == col of Pot */ 1204 apnz = api[row+1] - api[row]; 1205 Jptr = apj + api[row]; 1206 /* add non-zero cols of AP into the sorted linked list lnk */ 1207 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1208 } 1209 nnz = lnk[0]; 1210 1211 /* If free space is not available, double the total space in the list */ 1212 if (current_space->local_remaining<nnz) { 1213 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1214 nspacedouble++; 1215 } 1216 1217 /* Copy data into free space, and zero out denserows */ 1218 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1219 1220 current_space->array += nnz; 1221 current_space->local_used += nnz; 1222 current_space->local_remaining -= nnz; 1223 1224 coi[i+1] = coi[i] + nnz; 1225 } 1226 1227 ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 1228 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1229 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 1230 if (afill_tmp > afill) afill = afill_tmp; 1231 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 1232 1233 /* (3) send j-array (coj) of Co to other processors */ 1234 /*--------------------------------------------------*/ 1235 ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 1236 len_s = merge->len_s; 1237 merge->nsend = 0; 1238 1239 1240 /* determine row ownership */ 1241 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1242 merge->rowmap->n = pn; 1243 merge->rowmap->bs = 1; 1244 1245 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1246 owners = merge->rowmap->range; 1247 1248 /* determine the number of messages to send, their lengths */ 1249 ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 1250 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1251 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1252 1253 proc = 0; 1254 for (i=0; i<pon; i++) { 1255 while (prmap[i] >= owners[proc+1]) proc++; 1256 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 1257 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 1258 } 1259 1260 len = 0; /* max length of buf_si[], see (4) */ 1261 owners_co[0] = 0; 1262 for (proc=0; proc<size; proc++) { 1263 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1264 if (len_s[proc]) { 1265 merge->nsend++; 1266 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 1267 len += len_si[proc]; 1268 } 1269 } 1270 1271 /* determine the number and length of messages to receive for coi and coj */ 1272 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1273 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1274 1275 /* post the Irecv and Isend of coj */ 1276 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1277 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1278 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1279 for (proc=0, k=0; proc<size; proc++) { 1280 if (!len_s[proc]) continue; 1281 i = owners_co[proc]; 1282 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1283 k++; 1284 } 1285 1286 /* receives and sends of coj are complete */ 1287 for (i=0; i<merge->nrecv; i++) { 1288 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1289 } 1290 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1291 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1292 1293 /* (4) send and recv coi */ 1294 /*-----------------------*/ 1295 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1296 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1297 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1298 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1299 for (proc=0,k=0; proc<size; proc++) { 1300 if (!len_s[proc]) continue; 1301 /* form outgoing message for i-structure: 1302 buf_si[0]: nrows to be sent 1303 [1:nrows]: row index (global) 1304 [nrows+1:2*nrows+1]: i-structure index 1305 */ 1306 /*-------------------------------------------*/ 1307 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 1308 buf_si_i = buf_si + nrows+1; 1309 buf_si[0] = nrows; 1310 buf_si_i[0] = 0; 1311 nrows = 0; 1312 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1313 nzi = coi[i+1] - coi[i]; 1314 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1315 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1316 nrows++; 1317 } 1318 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1319 k++; 1320 buf_si += len_si[proc]; 1321 } 1322 i = merge->nrecv; 1323 while (i--) { 1324 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1325 } 1326 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1327 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1328 1329 ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 1330 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1331 ierr = PetscFree(swaits);CHKERRQ(ierr); 1332 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1333 1334 /* (5) compute the local portion of C (mpi mat) */ 1335 /*----------------------------------------------*/ 1336 ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 1337 1338 /* allocate pti array and free space for accumulating nonzero column info */ 1339 ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 1340 pti[0] = 0; 1341 1342 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1343 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am])); 1344 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1345 current_space = free_space; 1346 1347 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1348 for (k=0; k<merge->nrecv; k++) { 1349 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1350 nrows = *buf_ri_k[k]; 1351 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1352 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1353 } 1354 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 1355 1356 for (i=0; i<pn; i++) { 1357 /* add pdt[i,:]*AP into lnk */ 1358 pnz = pdti[i+1] - pdti[i]; 1359 ptJ = pdtj + pdti[i]; 1360 for (j=0; j<pnz; j++) { 1361 row = ptJ[j]; /* row of AP == col of Pt */ 1362 apnz = api[row+1] - api[row]; 1363 Jptr = apj + api[row]; 1364 /* add non-zero cols of AP into the sorted linked list lnk */ 1365 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1366 } 1367 1368 /* add received col data into lnk */ 1369 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1370 if (i == *nextrow[k]) { /* i-th row */ 1371 nzi = *(nextci[k]+1) - *nextci[k]; 1372 Jptr = buf_rj[k] + *nextci[k]; 1373 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1374 nextrow[k]++; nextci[k]++; 1375 } 1376 } 1377 nnz = lnk[0]; 1378 1379 /* if free space is not available, make more free space */ 1380 if (current_space->local_remaining<nnz) { 1381 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1382 nspacedouble++; 1383 } 1384 /* copy data into free space, then initialize lnk */ 1385 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1386 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1387 1388 current_space->array += nnz; 1389 current_space->local_used += nnz; 1390 current_space->local_remaining -= nnz; 1391 1392 pti[i+1] = pti[i] + nnz; 1393 } 1394 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 1395 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1396 1397 ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 1398 ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 1399 afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 1400 if (afill_tmp > afill) afill = afill_tmp; 1401 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1402 1403 /* (6) create symbolic parallel matrix Cmpi */ 1404 /*------------------------------------------*/ 1405 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1406 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1407 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 1408 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1409 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1410 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1411 1412 merge->bi = pti; /* Cseq->i */ 1413 merge->bj = ptj; /* Cseq->j */ 1414 merge->coi = coi; /* Co->i */ 1415 merge->coj = coj; /* Co->j */ 1416 merge->buf_ri = buf_ri; 1417 merge->buf_rj = buf_rj; 1418 merge->owners_co = owners_co; 1419 merge->destroy = Cmpi->ops->destroy; 1420 1421 /* attach the supporting struct to Cmpi for reuse */ 1422 c = (Mat_MPIAIJ*)Cmpi->data; 1423 c->ptap = ptap; 1424 ptap->api = api; 1425 ptap->apj = apj; 1426 ptap->duplicate = Cmpi->ops->duplicate; 1427 ptap->destroy = Cmpi->ops->destroy; 1428 1429 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1430 Cmpi->assembled = PETSC_FALSE; 1431 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1432 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1433 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap; 1434 *C = Cmpi; 1435 1436 /* flag 'scalable' determines which implementations to be used: 1437 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 1438 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 1439 /* set default scalable */ 1440 ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */ 1441 1442 ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 1443 if (!ptap->scalable) { /* Do dense axpy */ 1444 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 1445 } else { 1446 ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 1447 } 1448 1449 #if defined(PETSC_USE_INFO) 1450 if (pti[pn] != 0) { 1451 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1452 ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1453 } else { 1454 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1455 } 1456 #endif 1457 PetscFunctionReturn(0); 1458 } 1459 1460 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C) 1461 { 1462 PetscErrorCode ierr; 1463 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1464 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1465 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1466 Mat_SeqAIJ *p_loc,*p_oth; 1467 Mat_PtAPMPI *ptap; 1468 Mat_Merge_SeqsToMPI *merge; 1469 PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 1470 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 1471 PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1472 MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1473 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1474 MPI_Comm comm; 1475 PetscMPIInt size,rank,taga,*len_s; 1476 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1477 PetscInt **buf_ri,**buf_rj; 1478 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1479 MPI_Request *s_waits,*r_waits; 1480 MPI_Status *status; 1481 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1482 PetscInt *api,*apj,*coi,*coj; 1483 PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 1484 PetscBool scalable; 1485 #if defined(PTAP_PROFILE) 1486 PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 1487 #endif 1488 1489 PetscFunctionBegin; 1490 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1491 #if defined(PTAP_PROFILE) 1492 ierr = PetscTime(&t0);CHKERRQ(ierr); 1493 #endif 1494 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1495 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1496 1497 ptap = c->ptap; 1498 if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX"); 1499 merge = ptap->merge; 1500 apa = ptap->apa; 1501 scalable = ptap->scalable; 1502 1503 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 1504 /*-----------------------------------------------------*/ 1505 if (ptap->reuse == MAT_INITIAL_MATRIX) { 1506 /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1507 ptap->reuse = MAT_REUSE_MATRIX; 1508 } else { /* update numerical values of P_oth and P_loc */ 1509 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1510 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1511 } 1512 #if defined(PTAP_PROFILE) 1513 ierr = PetscTime(&t1);CHKERRQ(ierr); 1514 eP = t1-t0; 1515 #endif 1516 /* 1517 printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 1518 a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 1519 ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 1520 ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 1521 */ 1522 1523 /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1524 /*--------------------------------------------------------------*/ 1525 /* get data from symbolic products */ 1526 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1527 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1528 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 1529 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 1530 1531 coi = merge->coi; coj = merge->coj; 1532 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1533 1534 bi = merge->bi; bj = merge->bj; 1535 owners = merge->rowmap->range; 1536 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 1537 1538 api = ptap->api; apj = ptap->apj; 1539 1540 if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1541 ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 1542 #if defined(PTAP_PROFILE) 1543 ierr = PetscTime(&t1);CHKERRQ(ierr); 1544 #endif 1545 for (i=0; i<am; i++) { 1546 #if defined(PTAP_PROFILE) 1547 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1548 #endif 1549 /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1550 /*------------------------------------------------------------*/ 1551 apJ = apj + api[i]; 1552 1553 /* diagonal portion of A */ 1554 anz = adi[i+1] - adi[i]; 1555 adj = ad->j + adi[i]; 1556 ada = ad->a + adi[i]; 1557 for (j=0; j<anz; j++) { 1558 row = adj[j]; 1559 pnz = pi_loc[row+1] - pi_loc[row]; 1560 pj = pj_loc + pi_loc[row]; 1561 pa = pa_loc + pi_loc[row]; 1562 1563 /* perform dense axpy */ 1564 valtmp = ada[j]; 1565 for (k=0; k<pnz; k++) { 1566 apa[pj[k]] += valtmp*pa[k]; 1567 } 1568 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1569 } 1570 1571 /* off-diagonal portion of A */ 1572 anz = aoi[i+1] - aoi[i]; 1573 aoj = ao->j + aoi[i]; 1574 aoa = ao->a + aoi[i]; 1575 for (j=0; j<anz; j++) { 1576 row = aoj[j]; 1577 pnz = pi_oth[row+1] - pi_oth[row]; 1578 pj = pj_oth + pi_oth[row]; 1579 pa = pa_oth + pi_oth[row]; 1580 1581 /* perform dense axpy */ 1582 valtmp = aoa[j]; 1583 for (k=0; k<pnz; k++) { 1584 apa[pj[k]] += valtmp*pa[k]; 1585 } 1586 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1587 } 1588 #if defined(PTAP_PROFILE) 1589 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1590 et2_AP += t2_1 - t2_0; 1591 #endif 1592 1593 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1594 /*--------------------------------------------------------------*/ 1595 apnz = api[i+1] - api[i]; 1596 /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1597 pnz = po->i[i+1] - po->i[i]; 1598 poJ = po->j + po->i[i]; 1599 pA = po->a + po->i[i]; 1600 for (j=0; j<pnz; j++) { 1601 row = poJ[j]; 1602 cnz = coi[row+1] - coi[row]; 1603 cj = coj + coi[row]; 1604 ca = coa + coi[row]; 1605 /* perform dense axpy */ 1606 valtmp = pA[j]; 1607 for (k=0; k<cnz; k++) { 1608 ca[k] += valtmp*apa[cj[k]]; 1609 } 1610 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1611 } 1612 /* put the value into Cd (diagonal part) */ 1613 pnz = pd->i[i+1] - pd->i[i]; 1614 pdJ = pd->j + pd->i[i]; 1615 pA = pd->a + pd->i[i]; 1616 for (j=0; j<pnz; j++) { 1617 row = pdJ[j]; 1618 cnz = bi[row+1] - bi[row]; 1619 cj = bj + bi[row]; 1620 ca = ba + bi[row]; 1621 /* perform dense axpy */ 1622 valtmp = pA[j]; 1623 for (k=0; k<cnz; k++) { 1624 ca[k] += valtmp*apa[cj[k]]; 1625 } 1626 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1627 } 1628 1629 /* zero the current row of A*P */ 1630 for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1631 #if defined(PTAP_PROFILE) 1632 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 1633 ePtAP += t2_2 - t2_1; 1634 #endif 1635 } 1636 } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1637 ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 1638 /*-----------------------------------------------------------------------------------------*/ 1639 pA=pa_loc; 1640 for (i=0; i<am; i++) { 1641 #if defined(PTAP_PROFILE) 1642 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1643 #endif 1644 /* form i-th sparse row of A*P */ 1645 apnz = api[i+1] - api[i]; 1646 apJ = apj + api[i]; 1647 /* diagonal portion of A */ 1648 anz = adi[i+1] - adi[i]; 1649 adj = ad->j + adi[i]; 1650 ada = ad->a + adi[i]; 1651 for (j=0; j<anz; j++) { 1652 row = adj[j]; 1653 pnz = pi_loc[row+1] - pi_loc[row]; 1654 pj = pj_loc + pi_loc[row]; 1655 pa = pa_loc + pi_loc[row]; 1656 valtmp = ada[j]; 1657 nextp = 0; 1658 for (k=0; nextp<pnz; k++) { 1659 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1660 apa[k] += valtmp*pa[nextp++]; 1661 } 1662 } 1663 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1664 } 1665 /* off-diagonal portion of A */ 1666 anz = aoi[i+1] - aoi[i]; 1667 aoj = ao->j + aoi[i]; 1668 aoa = ao->a + aoi[i]; 1669 for (j=0; j<anz; j++) { 1670 row = aoj[j]; 1671 pnz = pi_oth[row+1] - pi_oth[row]; 1672 pj = pj_oth + pi_oth[row]; 1673 pa = pa_oth + pi_oth[row]; 1674 valtmp = aoa[j]; 1675 nextp = 0; 1676 for (k=0; nextp<pnz; k++) { 1677 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1678 apa[k] += valtmp*pa[nextp++]; 1679 } 1680 } 1681 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1682 } 1683 #if defined(PTAP_PROFILE) 1684 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1685 et2_AP += t2_1 - t2_0; 1686 #endif 1687 1688 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1689 /*--------------------------------------------------------------*/ 1690 pnz = pi_loc[i+1] - pi_loc[i]; 1691 pJ = pj_loc + pi_loc[i]; 1692 for (j=0; j<pnz; j++) { 1693 nextap = 0; 1694 row = pJ[j]; /* global index */ 1695 if (row < pcstart || row >=pcend) { /* put the value into Co */ 1696 row = *poJ; 1697 cj = coj + coi[row]; 1698 ca = coa + coi[row]; poJ++; 1699 } else { /* put the value into Cd */ 1700 row = *pdJ; 1701 cj = bj + bi[row]; 1702 ca = ba + bi[row]; pdJ++; 1703 } 1704 valtmp = pA[j]; 1705 for (k=0; nextap<apnz; k++) { 1706 if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 1707 } 1708 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1709 } 1710 pA += pnz; 1711 /* zero the current row info for A*P */ 1712 ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1713 #if defined(PTAP_PROFILE) 1714 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 1715 ePtAP += t2_2 - t2_1; 1716 #endif 1717 } 1718 } 1719 #if defined(PTAP_PROFILE) 1720 ierr = PetscTime(&t2);CHKERRQ(ierr); 1721 #endif 1722 1723 /* 3) send and recv matrix values coa */ 1724 /*------------------------------------*/ 1725 buf_ri = merge->buf_ri; 1726 buf_rj = merge->buf_rj; 1727 len_s = merge->len_s; 1728 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1729 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1730 1731 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1732 for (proc=0,k=0; proc<size; proc++) { 1733 if (!len_s[proc]) continue; 1734 i = merge->owners_co[proc]; 1735 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1736 k++; 1737 } 1738 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1739 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1740 1741 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1742 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1743 ierr = PetscFree(coa);CHKERRQ(ierr); 1744 #if defined(PTAP_PROFILE) 1745 ierr = PetscTime(&t3);CHKERRQ(ierr); 1746 #endif 1747 1748 /* 4) insert local Cseq and received values into Cmpi */ 1749 /*------------------------------------------------------*/ 1750 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1751 for (k=0; k<merge->nrecv; k++) { 1752 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1753 nrows = *(buf_ri_k[k]); 1754 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1755 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1756 } 1757 1758 for (i=0; i<cm; i++) { 1759 row = owners[rank] + i; /* global row index of C_seq */ 1760 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1761 ba_i = ba + bi[i]; 1762 bnz = bi[i+1] - bi[i]; 1763 /* add received vals into ba */ 1764 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1765 /* i-th row */ 1766 if (i == *nextrow[k]) { 1767 cnz = *(nextci[k]+1) - *nextci[k]; 1768 cj = buf_rj[k] + *(nextci[k]); 1769 ca = abuf_r[k] + *(nextci[k]); 1770 nextcj = 0; 1771 for (j=0; nextcj<cnz; j++) { 1772 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1773 ba_i[j] += ca[nextcj++]; 1774 } 1775 } 1776 nextrow[k]++; nextci[k]++; 1777 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1778 } 1779 } 1780 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1781 } 1782 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1783 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1784 1785 ierr = PetscFree(ba);CHKERRQ(ierr); 1786 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1787 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1788 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1789 #if defined(PTAP_PROFILE) 1790 ierr = PetscTime(&t4);CHKERRQ(ierr); 1791 if (rank==1) { 1792 ierr = PetscPrintf(MPI_COMM_SELF," [%d] PtAPNum %g/P + %g/PtAP( %g/A*P + %g/Pt*AP ) + %g/comm + %g/Cloc = %g\n\n",rank,eP,t2-t1,et2_AP,ePtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr); 1793 } 1794 #endif 1795 PetscFunctionReturn(0); 1796 } 1797