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