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,"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,"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 ierr = PetscInfo2(P,"Use scalable algorithm, pN %D, fill*nz_allocated %g\n",pN,fill*nz_local);CHKERRQ(ierr); 603 } 604 } 605 /* ierr = PetscInfo2(P,"Use algorithm %D, pN %D\n",alg,pN);CHKERRQ(ierr); */ 606 607 if (alg == 0) { 608 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr); 609 (*C)->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 610 PetscFunctionReturn(0); 611 612 #if defined(PETSC_HAVE_HYPRE) 613 } else if (alg == 2) { 614 /* Use boomerAMGBuildCoarseOperator */ 615 ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 616 PetscFunctionReturn(0); 617 #endif 618 } 619 620 #if defined(PTAP_PROFILE) 621 ierr = PetscTime(&t0);CHKERRQ(ierr); 622 #endif 623 624 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 625 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 626 627 if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 628 629 /* create symbolic parallel matrix Cmpi */ 630 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 631 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 632 633 /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 634 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 635 636 /* create struct Mat_PtAPMPI and attached it to C later */ 637 ierr = PetscNew(&ptap);CHKERRQ(ierr); 638 ptap->reuse = MAT_INITIAL_MATRIX; 639 640 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 641 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 642 /* get P_loc by taking all local rows of P */ 643 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 644 645 /* (0) compute Rd = Pd^T, Ro = Po^T */ 646 /* --------------------------------- */ 647 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 648 ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 649 #if defined(PTAP_PROFILE) 650 ierr = PetscTime(&t11);CHKERRQ(ierr); 651 #endif 652 653 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 654 /* ----------------------------------------------------------------- */ 655 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 656 if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 657 658 /* create and initialize a linked list */ 659 ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 660 MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 661 MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 662 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 663 /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */ 664 665 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 666 667 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 668 if (ao) { 669 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 670 } else { 671 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 672 } 673 current_space = free_space; 674 nspacedouble = 0; 675 676 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 677 api[0] = 0; 678 for (i=0; i<am; i++) { 679 /* diagonal portion: Ad[i,:]*P */ 680 ai = ad->i; pi = p_loc->i; 681 nzi = ai[i+1] - ai[i]; 682 aj = ad->j + ai[i]; 683 for (j=0; j<nzi; j++) { 684 row = aj[j]; 685 pnz = pi[row+1] - pi[row]; 686 Jptr = p_loc->j + pi[row]; 687 /* add non-zero cols of P into the sorted linked list lnk */ 688 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 689 } 690 /* off-diagonal portion: Ao[i,:]*P */ 691 if (ao) { 692 ai = ao->i; pi = p_oth->i; 693 nzi = ai[i+1] - ai[i]; 694 aj = ao->j + ai[i]; 695 for (j=0; j<nzi; j++) { 696 row = aj[j]; 697 pnz = pi[row+1] - pi[row]; 698 Jptr = p_oth->j + pi[row]; 699 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 700 } 701 } 702 apnz = lnk[0]; 703 api[i+1] = api[i] + apnz; 704 if (ap_rmax < apnz) ap_rmax = apnz; 705 706 /* if free space is not available, double the total space in the list */ 707 if (current_space->local_remaining<apnz) { 708 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 709 nspacedouble++; 710 } 711 712 /* Copy data into free space, then initialize lnk */ 713 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 714 715 current_space->array += apnz; 716 current_space->local_used += apnz; 717 current_space->local_remaining -= apnz; 718 } 719 /* Allocate space for apj and apv, initialize apj, and */ 720 /* destroy list of free space and other temporary array(s) */ 721 ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 722 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 723 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 724 725 /* Create AP_loc for reuse */ 726 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 727 728 #if defined(PETSC_USE_INFO) 729 if (ao) { 730 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 731 } else { 732 apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 733 } 734 ptap->AP_loc->info.mallocs = nspacedouble; 735 ptap->AP_loc->info.fill_ratio_given = fill; 736 ptap->AP_loc->info.fill_ratio_needed = apfill; 737 738 if (api[am]) { 739 ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 740 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 741 } else { 742 ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr); 743 } 744 #endif 745 746 #if defined(PTAP_PROFILE) 747 ierr = PetscTime(&t12);CHKERRQ(ierr); 748 #endif 749 750 /* (2-1) compute symbolic Co = Ro*AP_loc */ 751 /* ------------------------------------ */ 752 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 753 #if defined(PTAP_PROFILE) 754 ierr = PetscTime(&t1);CHKERRQ(ierr); 755 #endif 756 757 /* (3) send coj of C_oth to other processors */ 758 /* ------------------------------------------ */ 759 /* determine row ownership */ 760 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 761 rowmap->n = pn; 762 rowmap->bs = 1; 763 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 764 owners = rowmap->range; 765 766 /* determine the number of messages to send, their lengths */ 767 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 768 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 769 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 770 771 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 772 coi = c_oth->i; coj = c_oth->j; 773 con = ptap->C_oth->rmap->n; 774 proc = 0; 775 for (i=0; i<con; i++) { 776 while (prmap[i] >= owners[proc+1]) proc++; 777 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 778 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 779 } 780 781 len = 0; /* max length of buf_si[], see (4) */ 782 owners_co[0] = 0; 783 nsend = 0; 784 for (proc=0; proc<size; proc++) { 785 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 786 if (len_s[proc]) { 787 nsend++; 788 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 789 len += len_si[proc]; 790 } 791 } 792 793 /* determine the number and length of messages to receive for coi and coj */ 794 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 795 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 796 797 /* post the Irecv and Isend of coj */ 798 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 799 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 800 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 801 for (proc=0, k=0; proc<size; proc++) { 802 if (!len_s[proc]) continue; 803 i = owners_co[proc]; 804 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 805 k++; 806 } 807 808 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 809 /* ---------------------------------------- */ 810 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 811 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 812 813 /* receives coj are complete */ 814 for (i=0; i<nrecv; i++) { 815 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 816 } 817 ierr = PetscFree(rwaits);CHKERRQ(ierr); 818 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 819 820 /* add received column indices into ta to update Crmax */ 821 for (k=0; k<nrecv; k++) {/* k-th received message */ 822 Jptr = buf_rj[k]; 823 for (j=0; j<len_r[k]; j++) { 824 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 825 } 826 } 827 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 828 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 829 830 /* (4) send and recv coi */ 831 /*-----------------------*/ 832 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 833 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 834 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 835 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 836 for (proc=0,k=0; proc<size; proc++) { 837 if (!len_s[proc]) continue; 838 /* form outgoing message for i-structure: 839 buf_si[0]: nrows to be sent 840 [1:nrows]: row index (global) 841 [nrows+1:2*nrows+1]: i-structure index 842 */ 843 /*-------------------------------------------*/ 844 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 845 buf_si_i = buf_si + nrows+1; 846 buf_si[0] = nrows; 847 buf_si_i[0] = 0; 848 nrows = 0; 849 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 850 nzi = coi[i+1] - coi[i]; 851 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 852 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 853 nrows++; 854 } 855 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 856 k++; 857 buf_si += len_si[proc]; 858 } 859 for (i=0; i<nrecv; i++) { 860 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 861 } 862 ierr = PetscFree(rwaits);CHKERRQ(ierr); 863 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 864 865 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 866 ierr = PetscFree(len_ri);CHKERRQ(ierr); 867 ierr = PetscFree(swaits);CHKERRQ(ierr); 868 ierr = PetscFree(buf_s);CHKERRQ(ierr); 869 #if defined(PTAP_PROFILE) 870 ierr = PetscTime(&t2);CHKERRQ(ierr); 871 #endif 872 /* (5) compute the local portion of Cmpi */ 873 /* ------------------------------------------ */ 874 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 875 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 876 current_space = free_space; 877 878 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 879 for (k=0; k<nrecv; k++) { 880 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 881 nrows = *buf_ri_k[k]; 882 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 883 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 884 } 885 886 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 887 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 888 for (i=0; i<pn; i++) { 889 /* add C_loc into Cmpi */ 890 nzi = c_loc->i[i+1] - c_loc->i[i]; 891 Jptr = c_loc->j + c_loc->i[i]; 892 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 893 894 /* add received col data into lnk */ 895 for (k=0; k<nrecv; k++) { /* k-th received message */ 896 if (i == *nextrow[k]) { /* i-th row */ 897 nzi = *(nextci[k]+1) - *nextci[k]; 898 Jptr = buf_rj[k] + *nextci[k]; 899 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 900 nextrow[k]++; nextci[k]++; 901 } 902 } 903 nzi = lnk[0]; 904 905 /* copy data into free space, then initialize lnk */ 906 ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 907 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 908 } 909 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 910 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 911 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 912 #if defined(PTAP_PROFILE) 913 ierr = PetscTime(&t3);CHKERRQ(ierr); 914 #endif 915 916 /* local sizes and preallocation */ 917 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 918 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 919 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 920 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 921 922 /* members in merge */ 923 ierr = PetscFree(id_r);CHKERRQ(ierr); 924 ierr = PetscFree(len_r);CHKERRQ(ierr); 925 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 926 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 927 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 928 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 929 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 930 931 /* attach the supporting struct to Cmpi for reuse */ 932 c = (Mat_MPIAIJ*)Cmpi->data; 933 c->ptap = ptap; 934 ptap->duplicate = Cmpi->ops->duplicate; 935 ptap->destroy = Cmpi->ops->destroy; 936 937 if (alg == 1) { 938 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 939 } 940 941 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 942 Cmpi->assembled = PETSC_FALSE; 943 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 944 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 945 *C = Cmpi; 946 947 #if defined(PTAP_PROFILE) 948 ierr = PetscTime(&t4);CHKERRQ(ierr); 949 if (rank == 1) { 950 printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0); 951 } 952 #endif 953 PetscFunctionReturn(0); 954 } 955 956 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 957 { 958 PetscErrorCode ierr; 959 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 960 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 961 Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 962 Mat_PtAPMPI *ptap = c->ptap; 963 Mat AP_loc,C_loc,C_oth; 964 PetscInt i,rstart,rend,cm,ncols,row; 965 PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 966 PetscScalar *apa; 967 const PetscInt *cols; 968 const PetscScalar *vals; 969 #if defined(PTAP_PROFILE) 970 PetscMPIInt rank; 971 MPI_Comm comm; 972 PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 973 #endif 974 975 PetscFunctionBegin; 976 ierr = MatZeroEntries(C);CHKERRQ(ierr); 977 978 /* 1) get R = Pd^T,Ro = Po^T */ 979 #if defined(PTAP_PROFILE) 980 ierr = PetscTime(&t0);CHKERRQ(ierr); 981 #endif 982 if (ptap->reuse == MAT_REUSE_MATRIX) { 983 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 984 ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 985 } 986 #if defined(PTAP_PROFILE) 987 ierr = PetscTime(&t1);CHKERRQ(ierr); 988 eR = t1 - t0; 989 #endif 990 991 /* 2) get AP_loc */ 992 AP_loc = ptap->AP_loc; 993 ap = (Mat_SeqAIJ*)AP_loc->data; 994 995 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 996 /*-----------------------------------------------------*/ 997 if (ptap->reuse == MAT_REUSE_MATRIX) { 998 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 999 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1000 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1001 } 1002 1003 /* 2-2) compute numeric A_loc*P - dominating part */ 1004 /* ---------------------------------------------- */ 1005 /* get data from symbolic products */ 1006 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1007 if (ptap->P_oth) { 1008 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1009 } 1010 apa = ptap->apa; 1011 api = ap->i; 1012 apj = ap->j; 1013 for (i=0; i<am; i++) { 1014 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 1015 AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 1016 apnz = api[i+1] - api[i]; 1017 for (j=0; j<apnz; j++) { 1018 col = apj[j+api[i]]; 1019 ap->a[j+ap->i[i]] = apa[col]; 1020 apa[col] = 0.0; 1021 } 1022 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1023 } 1024 #if defined(PTAP_PROFILE) 1025 ierr = PetscTime(&t2);CHKERRQ(ierr); 1026 eAP = t2 - t1; 1027 #endif 1028 1029 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 1030 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 1031 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 1032 C_loc = ptap->C_loc; 1033 C_oth = ptap->C_oth; 1034 #if defined(PTAP_PROFILE) 1035 ierr = PetscTime(&t3);CHKERRQ(ierr); 1036 eCseq = t3 - t2; 1037 #endif 1038 1039 /* add C_loc and Co to to C */ 1040 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 1041 1042 /* C_loc -> C */ 1043 cm = C_loc->rmap->N; 1044 c_seq = (Mat_SeqAIJ*)C_loc->data; 1045 cols = c_seq->j; 1046 vals = c_seq->a; 1047 for (i=0; i<cm; i++) { 1048 ncols = c_seq->i[i+1] - c_seq->i[i]; 1049 row = rstart + i; 1050 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1051 cols += ncols; vals += ncols; 1052 } 1053 1054 /* Co -> C, off-processor part */ 1055 cm = C_oth->rmap->N; 1056 c_seq = (Mat_SeqAIJ*)C_oth->data; 1057 cols = c_seq->j; 1058 vals = c_seq->a; 1059 for (i=0; i<cm; i++) { 1060 ncols = c_seq->i[i+1] - c_seq->i[i]; 1061 row = p->garray[i]; 1062 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1063 cols += ncols; vals += ncols; 1064 } 1065 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1066 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1067 #if defined(PTAP_PROFILE) 1068 ierr = PetscTime(&t4);CHKERRQ(ierr); 1069 eCmpi = t4 - t3; 1070 1071 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1072 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1073 if (rank==1) { 1074 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); 1075 } 1076 #endif 1077 ptap->reuse = MAT_REUSE_MATRIX; 1078 PetscFunctionReturn(0); 1079 } 1080 1081 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C) 1082 { 1083 PetscErrorCode ierr; 1084 Mat Cmpi; 1085 Mat_PtAPMPI *ptap; 1086 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1087 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 1088 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1089 Mat_SeqAIJ *p_loc,*p_oth; 1090 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 1091 PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 1092 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1093 PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 1094 PetscBT lnkbt; 1095 MPI_Comm comm; 1096 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 1097 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1098 PetscInt len,proc,*dnz,*onz,*owners; 1099 PetscInt nzi,*pti,*ptj; 1100 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1101 MPI_Request *swaits,*rwaits; 1102 MPI_Status *sstatus,rstatus; 1103 Mat_Merge_SeqsToMPI *merge; 1104 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 1105 PetscReal afill=1.0,afill_tmp; 1106 1107 PetscFunctionBegin; 1108 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1109 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1110 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1111 1112 /* create struct Mat_PtAPMPI and attached it to C later */ 1113 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1114 ierr = PetscNew(&merge);CHKERRQ(ierr); 1115 ptap->merge = merge; 1116 ptap->reuse = MAT_INITIAL_MATRIX; 1117 1118 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1119 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1120 1121 /* get P_loc by taking all local rows of P */ 1122 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1123 1124 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1125 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1126 pi_loc = p_loc->i; pj_loc = p_loc->j; 1127 pi_oth = p_oth->i; pj_oth = p_oth->j; 1128 1129 /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 1130 /*--------------------------------------------------------------------------*/ 1131 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 1132 api[0] = 0; 1133 1134 /* create and initialize a linked list */ 1135 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1136 1137 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 1138 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 1139 current_space = free_space; 1140 1141 for (i=0; i<am; i++) { 1142 /* diagonal portion of A */ 1143 nzi = adi[i+1] - adi[i]; 1144 aj = ad->j + adi[i]; 1145 for (j=0; j<nzi; j++) { 1146 row = aj[j]; 1147 pnz = pi_loc[row+1] - pi_loc[row]; 1148 Jptr = pj_loc + pi_loc[row]; 1149 /* add non-zero cols of P into the sorted linked list lnk */ 1150 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1151 } 1152 /* off-diagonal portion of A */ 1153 nzi = aoi[i+1] - aoi[i]; 1154 aj = ao->j + aoi[i]; 1155 for (j=0; j<nzi; j++) { 1156 row = aj[j]; 1157 pnz = pi_oth[row+1] - pi_oth[row]; 1158 Jptr = pj_oth + pi_oth[row]; 1159 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1160 } 1161 apnz = lnk[0]; 1162 api[i+1] = api[i] + apnz; 1163 if (ap_rmax < apnz) ap_rmax = apnz; 1164 1165 /* if free space is not available, double the total space in the list */ 1166 if (current_space->local_remaining<apnz) { 1167 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1168 nspacedouble++; 1169 } 1170 1171 /* Copy data into free space, then initialize lnk */ 1172 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1173 1174 current_space->array += apnz; 1175 current_space->local_used += apnz; 1176 current_space->local_remaining -= apnz; 1177 } 1178 1179 /* Allocate space for apj, initialize apj, and */ 1180 /* destroy list of free space and other temporary array(s) */ 1181 ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 1182 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 1183 afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 1184 if (afill_tmp > afill) afill = afill_tmp; 1185 1186 /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 1187 /*-----------------------------------------------------------------*/ 1188 ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 1189 1190 /* then, compute symbolic Co = (p->B)^T*AP */ 1191 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1192 >= (num of nonzero rows of C_seq) - pn */ 1193 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1194 coi[0] = 0; 1195 1196 /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 1197 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am])); 1198 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1199 current_space = free_space; 1200 1201 for (i=0; i<pon; i++) { 1202 pnz = poti[i+1] - poti[i]; 1203 ptJ = potj + poti[i]; 1204 for (j=0; j<pnz; j++) { 1205 row = ptJ[j]; /* row of AP == col of Pot */ 1206 apnz = api[row+1] - api[row]; 1207 Jptr = apj + api[row]; 1208 /* add non-zero cols of AP into the sorted linked list lnk */ 1209 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1210 } 1211 nnz = lnk[0]; 1212 1213 /* If free space is not available, double the total space in the list */ 1214 if (current_space->local_remaining<nnz) { 1215 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1216 nspacedouble++; 1217 } 1218 1219 /* Copy data into free space, and zero out denserows */ 1220 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1221 1222 current_space->array += nnz; 1223 current_space->local_used += nnz; 1224 current_space->local_remaining -= nnz; 1225 1226 coi[i+1] = coi[i] + nnz; 1227 } 1228 1229 ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 1230 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1231 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 1232 if (afill_tmp > afill) afill = afill_tmp; 1233 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 1234 1235 /* (3) send j-array (coj) of Co to other processors */ 1236 /*--------------------------------------------------*/ 1237 ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 1238 len_s = merge->len_s; 1239 merge->nsend = 0; 1240 1241 1242 /* determine row ownership */ 1243 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1244 merge->rowmap->n = pn; 1245 merge->rowmap->bs = 1; 1246 1247 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1248 owners = merge->rowmap->range; 1249 1250 /* determine the number of messages to send, their lengths */ 1251 ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 1252 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1253 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1254 1255 proc = 0; 1256 for (i=0; i<pon; i++) { 1257 while (prmap[i] >= owners[proc+1]) proc++; 1258 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 1259 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 1260 } 1261 1262 len = 0; /* max length of buf_si[], see (4) */ 1263 owners_co[0] = 0; 1264 for (proc=0; proc<size; proc++) { 1265 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1266 if (len_s[proc]) { 1267 merge->nsend++; 1268 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 1269 len += len_si[proc]; 1270 } 1271 } 1272 1273 /* determine the number and length of messages to receive for coi and coj */ 1274 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1275 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1276 1277 /* post the Irecv and Isend of coj */ 1278 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1279 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1280 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1281 for (proc=0, k=0; proc<size; proc++) { 1282 if (!len_s[proc]) continue; 1283 i = owners_co[proc]; 1284 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1285 k++; 1286 } 1287 1288 /* receives and sends of coj are complete */ 1289 for (i=0; i<merge->nrecv; i++) { 1290 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1291 } 1292 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1293 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1294 1295 /* (4) send and recv coi */ 1296 /*-----------------------*/ 1297 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1298 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1299 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1300 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1301 for (proc=0,k=0; proc<size; proc++) { 1302 if (!len_s[proc]) continue; 1303 /* form outgoing message for i-structure: 1304 buf_si[0]: nrows to be sent 1305 [1:nrows]: row index (global) 1306 [nrows+1:2*nrows+1]: i-structure index 1307 */ 1308 /*-------------------------------------------*/ 1309 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 1310 buf_si_i = buf_si + nrows+1; 1311 buf_si[0] = nrows; 1312 buf_si_i[0] = 0; 1313 nrows = 0; 1314 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1315 nzi = coi[i+1] - coi[i]; 1316 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1317 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1318 nrows++; 1319 } 1320 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1321 k++; 1322 buf_si += len_si[proc]; 1323 } 1324 i = merge->nrecv; 1325 while (i--) { 1326 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1327 } 1328 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1329 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1330 1331 ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 1332 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1333 ierr = PetscFree(swaits);CHKERRQ(ierr); 1334 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1335 1336 /* (5) compute the local portion of C (mpi mat) */ 1337 /*----------------------------------------------*/ 1338 ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 1339 1340 /* allocate pti array and free space for accumulating nonzero column info */ 1341 ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 1342 pti[0] = 0; 1343 1344 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1345 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am])); 1346 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1347 current_space = free_space; 1348 1349 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1350 for (k=0; k<merge->nrecv; k++) { 1351 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1352 nrows = *buf_ri_k[k]; 1353 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1354 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1355 } 1356 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 1357 1358 for (i=0; i<pn; i++) { 1359 /* add pdt[i,:]*AP into lnk */ 1360 pnz = pdti[i+1] - pdti[i]; 1361 ptJ = pdtj + pdti[i]; 1362 for (j=0; j<pnz; j++) { 1363 row = ptJ[j]; /* row of AP == col of Pt */ 1364 apnz = api[row+1] - api[row]; 1365 Jptr = apj + api[row]; 1366 /* add non-zero cols of AP into the sorted linked list lnk */ 1367 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1368 } 1369 1370 /* add received col data into lnk */ 1371 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1372 if (i == *nextrow[k]) { /* i-th row */ 1373 nzi = *(nextci[k]+1) - *nextci[k]; 1374 Jptr = buf_rj[k] + *nextci[k]; 1375 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1376 nextrow[k]++; nextci[k]++; 1377 } 1378 } 1379 nnz = lnk[0]; 1380 1381 /* if free space is not available, make more free space */ 1382 if (current_space->local_remaining<nnz) { 1383 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1384 nspacedouble++; 1385 } 1386 /* copy data into free space, then initialize lnk */ 1387 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1388 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1389 1390 current_space->array += nnz; 1391 current_space->local_used += nnz; 1392 current_space->local_remaining -= nnz; 1393 1394 pti[i+1] = pti[i] + nnz; 1395 } 1396 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 1397 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1398 1399 ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 1400 ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 1401 afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 1402 if (afill_tmp > afill) afill = afill_tmp; 1403 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1404 1405 /* (6) create symbolic parallel matrix Cmpi */ 1406 /*------------------------------------------*/ 1407 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1408 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1409 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 1410 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1411 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1412 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1413 1414 merge->bi = pti; /* Cseq->i */ 1415 merge->bj = ptj; /* Cseq->j */ 1416 merge->coi = coi; /* Co->i */ 1417 merge->coj = coj; /* Co->j */ 1418 merge->buf_ri = buf_ri; 1419 merge->buf_rj = buf_rj; 1420 merge->owners_co = owners_co; 1421 merge->destroy = Cmpi->ops->destroy; 1422 1423 /* attach the supporting struct to Cmpi for reuse */ 1424 c = (Mat_MPIAIJ*)Cmpi->data; 1425 c->ptap = ptap; 1426 ptap->api = api; 1427 ptap->apj = apj; 1428 ptap->duplicate = Cmpi->ops->duplicate; 1429 ptap->destroy = Cmpi->ops->destroy; 1430 1431 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1432 Cmpi->assembled = PETSC_FALSE; 1433 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1434 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1435 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap; 1436 *C = Cmpi; 1437 1438 /* flag 'scalable' determines which implementations to be used: 1439 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 1440 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 1441 /* set default scalable */ 1442 ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */ 1443 1444 ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 1445 if (!ptap->scalable) { /* Do dense axpy */ 1446 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 1447 } else { 1448 ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 1449 } 1450 1451 #if defined(PETSC_USE_INFO) 1452 if (pti[pn] != 0) { 1453 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1454 ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1455 } else { 1456 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1457 } 1458 #endif 1459 PetscFunctionReturn(0); 1460 } 1461 1462 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C) 1463 { 1464 PetscErrorCode ierr; 1465 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1466 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1467 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1468 Mat_SeqAIJ *p_loc,*p_oth; 1469 Mat_PtAPMPI *ptap; 1470 Mat_Merge_SeqsToMPI *merge; 1471 PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 1472 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 1473 PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1474 MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1475 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1476 MPI_Comm comm; 1477 PetscMPIInt size,rank,taga,*len_s; 1478 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1479 PetscInt **buf_ri,**buf_rj; 1480 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1481 MPI_Request *s_waits,*r_waits; 1482 MPI_Status *status; 1483 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1484 PetscInt *api,*apj,*coi,*coj; 1485 PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 1486 PetscBool scalable; 1487 #if defined(PTAP_PROFILE) 1488 PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 1489 #endif 1490 1491 PetscFunctionBegin; 1492 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1493 #if defined(PTAP_PROFILE) 1494 ierr = PetscTime(&t0);CHKERRQ(ierr); 1495 #endif 1496 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1497 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1498 1499 ptap = c->ptap; 1500 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"); 1501 merge = ptap->merge; 1502 apa = ptap->apa; 1503 scalable = ptap->scalable; 1504 1505 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 1506 /*-----------------------------------------------------*/ 1507 if (ptap->reuse == MAT_INITIAL_MATRIX) { 1508 /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1509 ptap->reuse = MAT_REUSE_MATRIX; 1510 } else { /* update numerical values of P_oth and P_loc */ 1511 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1512 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1513 } 1514 #if defined(PTAP_PROFILE) 1515 ierr = PetscTime(&t1);CHKERRQ(ierr); 1516 eP = t1-t0; 1517 #endif 1518 /* 1519 printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 1520 a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 1521 ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 1522 ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 1523 */ 1524 1525 /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1526 /*--------------------------------------------------------------*/ 1527 /* get data from symbolic products */ 1528 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1529 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1530 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 1531 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 1532 1533 coi = merge->coi; coj = merge->coj; 1534 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1535 1536 bi = merge->bi; bj = merge->bj; 1537 owners = merge->rowmap->range; 1538 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 1539 1540 api = ptap->api; apj = ptap->apj; 1541 1542 if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1543 ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 1544 #if defined(PTAP_PROFILE) 1545 ierr = PetscTime(&t1);CHKERRQ(ierr); 1546 #endif 1547 for (i=0; i<am; i++) { 1548 #if defined(PTAP_PROFILE) 1549 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1550 #endif 1551 /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1552 /*------------------------------------------------------------*/ 1553 apJ = apj + api[i]; 1554 1555 /* diagonal portion of A */ 1556 anz = adi[i+1] - adi[i]; 1557 adj = ad->j + adi[i]; 1558 ada = ad->a + adi[i]; 1559 for (j=0; j<anz; j++) { 1560 row = adj[j]; 1561 pnz = pi_loc[row+1] - pi_loc[row]; 1562 pj = pj_loc + pi_loc[row]; 1563 pa = pa_loc + pi_loc[row]; 1564 1565 /* perform dense axpy */ 1566 valtmp = ada[j]; 1567 for (k=0; k<pnz; k++) { 1568 apa[pj[k]] += valtmp*pa[k]; 1569 } 1570 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1571 } 1572 1573 /* off-diagonal portion of A */ 1574 anz = aoi[i+1] - aoi[i]; 1575 aoj = ao->j + aoi[i]; 1576 aoa = ao->a + aoi[i]; 1577 for (j=0; j<anz; j++) { 1578 row = aoj[j]; 1579 pnz = pi_oth[row+1] - pi_oth[row]; 1580 pj = pj_oth + pi_oth[row]; 1581 pa = pa_oth + pi_oth[row]; 1582 1583 /* perform dense axpy */ 1584 valtmp = aoa[j]; 1585 for (k=0; k<pnz; k++) { 1586 apa[pj[k]] += valtmp*pa[k]; 1587 } 1588 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1589 } 1590 #if defined(PTAP_PROFILE) 1591 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1592 et2_AP += t2_1 - t2_0; 1593 #endif 1594 1595 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1596 /*--------------------------------------------------------------*/ 1597 apnz = api[i+1] - api[i]; 1598 /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1599 pnz = po->i[i+1] - po->i[i]; 1600 poJ = po->j + po->i[i]; 1601 pA = po->a + po->i[i]; 1602 for (j=0; j<pnz; j++) { 1603 row = poJ[j]; 1604 cnz = coi[row+1] - coi[row]; 1605 cj = coj + coi[row]; 1606 ca = coa + coi[row]; 1607 /* perform dense axpy */ 1608 valtmp = pA[j]; 1609 for (k=0; k<cnz; k++) { 1610 ca[k] += valtmp*apa[cj[k]]; 1611 } 1612 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1613 } 1614 /* put the value into Cd (diagonal part) */ 1615 pnz = pd->i[i+1] - pd->i[i]; 1616 pdJ = pd->j + pd->i[i]; 1617 pA = pd->a + pd->i[i]; 1618 for (j=0; j<pnz; j++) { 1619 row = pdJ[j]; 1620 cnz = bi[row+1] - bi[row]; 1621 cj = bj + bi[row]; 1622 ca = ba + bi[row]; 1623 /* perform dense axpy */ 1624 valtmp = pA[j]; 1625 for (k=0; k<cnz; k++) { 1626 ca[k] += valtmp*apa[cj[k]]; 1627 } 1628 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1629 } 1630 1631 /* zero the current row of A*P */ 1632 for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1633 #if defined(PTAP_PROFILE) 1634 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 1635 ePtAP += t2_2 - t2_1; 1636 #endif 1637 } 1638 } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1639 ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 1640 /*-----------------------------------------------------------------------------------------*/ 1641 pA=pa_loc; 1642 for (i=0; i<am; i++) { 1643 #if defined(PTAP_PROFILE) 1644 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1645 #endif 1646 /* form i-th sparse row of A*P */ 1647 apnz = api[i+1] - api[i]; 1648 apJ = apj + api[i]; 1649 /* diagonal portion of A */ 1650 anz = adi[i+1] - adi[i]; 1651 adj = ad->j + adi[i]; 1652 ada = ad->a + adi[i]; 1653 for (j=0; j<anz; j++) { 1654 row = adj[j]; 1655 pnz = pi_loc[row+1] - pi_loc[row]; 1656 pj = pj_loc + pi_loc[row]; 1657 pa = pa_loc + pi_loc[row]; 1658 valtmp = ada[j]; 1659 nextp = 0; 1660 for (k=0; nextp<pnz; k++) { 1661 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1662 apa[k] += valtmp*pa[nextp++]; 1663 } 1664 } 1665 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1666 } 1667 /* off-diagonal portion of A */ 1668 anz = aoi[i+1] - aoi[i]; 1669 aoj = ao->j + aoi[i]; 1670 aoa = ao->a + aoi[i]; 1671 for (j=0; j<anz; j++) { 1672 row = aoj[j]; 1673 pnz = pi_oth[row+1] - pi_oth[row]; 1674 pj = pj_oth + pi_oth[row]; 1675 pa = pa_oth + pi_oth[row]; 1676 valtmp = aoa[j]; 1677 nextp = 0; 1678 for (k=0; nextp<pnz; k++) { 1679 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1680 apa[k] += valtmp*pa[nextp++]; 1681 } 1682 } 1683 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1684 } 1685 #if defined(PTAP_PROFILE) 1686 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1687 et2_AP += t2_1 - t2_0; 1688 #endif 1689 1690 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1691 /*--------------------------------------------------------------*/ 1692 pnz = pi_loc[i+1] - pi_loc[i]; 1693 pJ = pj_loc + pi_loc[i]; 1694 for (j=0; j<pnz; j++) { 1695 nextap = 0; 1696 row = pJ[j]; /* global index */ 1697 if (row < pcstart || row >=pcend) { /* put the value into Co */ 1698 row = *poJ; 1699 cj = coj + coi[row]; 1700 ca = coa + coi[row]; poJ++; 1701 } else { /* put the value into Cd */ 1702 row = *pdJ; 1703 cj = bj + bi[row]; 1704 ca = ba + bi[row]; pdJ++; 1705 } 1706 valtmp = pA[j]; 1707 for (k=0; nextap<apnz; k++) { 1708 if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 1709 } 1710 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1711 } 1712 pA += pnz; 1713 /* zero the current row info for A*P */ 1714 ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1715 #if defined(PTAP_PROFILE) 1716 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 1717 ePtAP += t2_2 - t2_1; 1718 #endif 1719 } 1720 } 1721 #if defined(PTAP_PROFILE) 1722 ierr = PetscTime(&t2);CHKERRQ(ierr); 1723 #endif 1724 1725 /* 3) send and recv matrix values coa */ 1726 /*------------------------------------*/ 1727 buf_ri = merge->buf_ri; 1728 buf_rj = merge->buf_rj; 1729 len_s = merge->len_s; 1730 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1731 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1732 1733 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1734 for (proc=0,k=0; proc<size; proc++) { 1735 if (!len_s[proc]) continue; 1736 i = merge->owners_co[proc]; 1737 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1738 k++; 1739 } 1740 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1741 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1742 1743 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1744 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1745 ierr = PetscFree(coa);CHKERRQ(ierr); 1746 #if defined(PTAP_PROFILE) 1747 ierr = PetscTime(&t3);CHKERRQ(ierr); 1748 #endif 1749 1750 /* 4) insert local Cseq and received values into Cmpi */ 1751 /*------------------------------------------------------*/ 1752 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1753 for (k=0; k<merge->nrecv; k++) { 1754 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1755 nrows = *(buf_ri_k[k]); 1756 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1757 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1758 } 1759 1760 for (i=0; i<cm; i++) { 1761 row = owners[rank] + i; /* global row index of C_seq */ 1762 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1763 ba_i = ba + bi[i]; 1764 bnz = bi[i+1] - bi[i]; 1765 /* add received vals into ba */ 1766 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1767 /* i-th row */ 1768 if (i == *nextrow[k]) { 1769 cnz = *(nextci[k]+1) - *nextci[k]; 1770 cj = buf_rj[k] + *(nextci[k]); 1771 ca = abuf_r[k] + *(nextci[k]); 1772 nextcj = 0; 1773 for (j=0; nextcj<cnz; j++) { 1774 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1775 ba_i[j] += ca[nextcj++]; 1776 } 1777 } 1778 nextrow[k]++; nextci[k]++; 1779 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1780 } 1781 } 1782 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1783 } 1784 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1785 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1786 1787 ierr = PetscFree(ba);CHKERRQ(ierr); 1788 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1789 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1790 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1791 #if defined(PTAP_PROFILE) 1792 ierr = PetscTime(&t4);CHKERRQ(ierr); 1793 if (rank==1) { 1794 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); 1795 } 1796 #endif 1797 PetscFunctionReturn(0); 1798 } 1799