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