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 #undef __FUNCT__ 17 #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP" 18 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 19 { 20 PetscErrorCode ierr; 21 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 22 Mat_PtAPMPI *ptap=a->ptap; 23 24 PetscFunctionBegin; 25 if (ptap) { 26 Mat_Merge_SeqsToMPI *merge=ptap->merge; 27 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 28 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 29 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 30 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 31 ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 32 ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 33 ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 34 if (ptap->AP_loc) { /* used by alg_rap */ 35 Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 36 ierr = PetscFree(ap->i);CHKERRQ(ierr); 37 ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 38 ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 39 } else { /* used by alg_ptap */ 40 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 41 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 42 } 43 ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 44 ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 45 if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 46 47 if (merge) { /* used by alg_ptap */ 48 ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 49 ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 50 ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 51 ierr = PetscFree(merge->bi);CHKERRQ(ierr); 52 ierr = PetscFree(merge->bj);CHKERRQ(ierr); 53 ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 54 ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 55 ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 56 ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 57 ierr = PetscFree(merge->coi);CHKERRQ(ierr); 58 ierr = PetscFree(merge->coj);CHKERRQ(ierr); 59 ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 60 ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 61 ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 62 } 63 64 ierr = ptap->destroy(A);CHKERRQ(ierr); 65 ierr = PetscFree(ptap);CHKERRQ(ierr); 66 } 67 PetscFunctionReturn(0); 68 } 69 70 #undef __FUNCT__ 71 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP" 72 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 73 { 74 PetscErrorCode ierr; 75 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 76 Mat_PtAPMPI *ptap = a->ptap; 77 78 PetscFunctionBegin; 79 ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 80 (*M)->ops->destroy = ptap->destroy; 81 (*M)->ops->duplicate = ptap->duplicate; 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 87 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 88 { 89 PetscErrorCode ierr; 90 PetscBool rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */ 91 MPI_Comm comm; 92 93 PetscFunctionBegin; 94 /* check if matrix local sizes are compatible */ 95 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 96 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); 97 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); 98 99 ierr = PetscOptionsGetBool(NULL,NULL,"-matrap",&rap,NULL);CHKERRQ(ierr); 100 if (scall == MAT_INITIAL_MATRIX) { 101 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 102 if (rap) { /* do R=P^T locally, then C=R*A*P */ 103 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 104 } else { /* do P^T*A*P */ 105 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr); 106 } 107 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 108 } 109 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 110 ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 111 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 112 PetscFunctionReturn(0); 113 } 114 115 #undef __FUNCT__ 116 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 117 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 118 { 119 PetscErrorCode ierr; 120 Mat_PtAPMPI *ptap; 121 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 122 MPI_Comm comm; 123 PetscMPIInt size,rank; 124 Mat Cmpi; 125 PetscFreeSpaceList free_space=NULL,current_space=NULL; 126 PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 127 PetscInt *lnk,i,k,pnz,row,nsend; 128 PetscBT lnkbt; 129 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 130 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 131 PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 132 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 133 MPI_Request *swaits,*rwaits; 134 MPI_Status *sstatus,rstatus; 135 PetscLayout rowmap; 136 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 137 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 138 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 139 Mat_SeqAIJ *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth; 140 PetscScalar *apv; 141 PetscTable ta; 142 const char *algTypes[2] = {"scalable","nonscalable"}; 143 PetscInt alg=1; /* set default algorithm */ 144 #if defined(PETSC_USE_INFO) 145 PetscReal apfill; 146 #endif 147 #if defined(PTAP_PROFILE) 148 PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 149 #endif 150 151 PetscFunctionBegin; 152 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 153 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 154 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 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,&ptap->P_oth);CHKERRQ(ierr); 165 /* get P_loc by taking all local rows of P */ 166 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 167 168 /* (0) compute Rd = Pd^T, Ro = Po^T */ 169 /* --------------------------------- */ 170 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 171 ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 172 #if defined(PTAP_PROFILE) 173 ierr = PetscTime(&t11);CHKERRQ(ierr); 174 #endif 175 176 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 177 /* ----------------------------------------------------------------- */ 178 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 179 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 180 181 /* create and initialize a linked list */ 182 ierr = PetscTableCreate(pN,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 183 MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 184 MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 185 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 186 /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */ 187 188 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 189 190 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 191 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 192 current_space = free_space; 193 nspacedouble = 0; 194 195 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 196 api[0] = 0; 197 for (i=0; i<am; i++) { 198 /* diagonal portion: Ad[i,:]*P */ 199 ai = ad->i; pi = p_loc->i; 200 nzi = ai[i+1] - ai[i]; 201 aj = ad->j + ai[i]; 202 for (j=0; j<nzi; j++) { 203 row = aj[j]; 204 pnz = pi[row+1] - pi[row]; 205 Jptr = p_loc->j + pi[row]; 206 /* add non-zero cols of P into the sorted linked list lnk */ 207 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 208 } 209 /* off-diagonal portion: Ao[i,:]*P */ 210 ai = ao->i; pi = p_oth->i; 211 nzi = ai[i+1] - ai[i]; 212 aj = ao->j + ai[i]; 213 for (j=0; j<nzi; j++) { 214 row = aj[j]; 215 pnz = pi[row+1] - pi[row]; 216 Jptr = p_oth->j + pi[row]; 217 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 218 } 219 apnz = lnk[0]; 220 api[i+1] = api[i] + apnz; 221 if (ap_rmax < apnz) ap_rmax = apnz; 222 223 /* if free space is not available, double the total space in the list */ 224 if (current_space->local_remaining<apnz) { 225 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 226 nspacedouble++; 227 } 228 229 /* Copy data into free space, then initialize lnk */ 230 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 231 232 current_space->array += apnz; 233 current_space->local_used += apnz; 234 current_space->local_remaining -= apnz; 235 } 236 /* Allocate space for apj and apv, initialize apj, and */ 237 /* destroy list of free space and other temporary array(s) */ 238 ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 239 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 240 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 241 242 /* Create AP_loc for reuse */ 243 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 244 245 #if defined(PETSC_USE_INFO) 246 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 247 ptap->AP_loc->info.mallocs = nspacedouble; 248 ptap->AP_loc->info.fill_ratio_given = fill; 249 ptap->AP_loc->info.fill_ratio_needed = apfill; 250 251 if (api[am]) { 252 ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 253 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 254 } else { 255 ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr); 256 } 257 #endif 258 259 #if defined(PTAP_PROFILE) 260 ierr = PetscTime(&t12);CHKERRQ(ierr); 261 #endif 262 263 /* (2-1) compute symbolic Co = Ro*AP_loc */ 264 /* ------------------------------------ */ 265 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 266 #if defined(PTAP_PROFILE) 267 ierr = PetscTime(&t1);CHKERRQ(ierr); 268 #endif 269 270 /* (3) send coj of C_oth to other processors */ 271 /* ------------------------------------------ */ 272 /* determine row ownership */ 273 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 274 rowmap->n = pn; 275 rowmap->bs = 1; 276 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 277 owners = rowmap->range; 278 279 /* determine the number of messages to send, their lengths */ 280 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 281 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 282 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 283 284 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 285 coi = c_oth->i; coj = c_oth->j; 286 con = ptap->C_oth->rmap->n; 287 proc = 0; 288 for (i=0; i<con; i++) { 289 while (prmap[i] >= owners[proc+1]) proc++; 290 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 291 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 292 } 293 294 len = 0; /* max length of buf_si[], see (4) */ 295 owners_co[0] = 0; 296 nsend = 0; 297 for (proc=0; proc<size; proc++) { 298 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 299 if (len_s[proc]) { 300 nsend++; 301 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 302 len += len_si[proc]; 303 } 304 } 305 306 /* determine the number and length of messages to receive for coi and coj */ 307 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 308 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 309 310 /* post the Irecv and Isend of coj */ 311 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 312 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 313 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 314 for (proc=0, k=0; proc<size; proc++) { 315 if (!len_s[proc]) continue; 316 i = owners_co[proc]; 317 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 318 k++; 319 } 320 321 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 322 /* ---------------------------------------- */ 323 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 324 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 325 326 /* receives coj are complete */ 327 for (i=0; i<nrecv; i++) { 328 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 329 } 330 ierr = PetscFree(rwaits);CHKERRQ(ierr); 331 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 332 333 /* add received column indices into ta to update Crmax */ 334 for (k=0; k<nrecv; k++) {/* k-th received message */ 335 Jptr = buf_rj[k]; 336 for (j=0; j<len_r[k]; j++) { 337 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 338 } 339 } 340 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 341 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 342 343 /* (4) send and recv coi */ 344 /*-----------------------*/ 345 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 346 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 347 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 348 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 349 for (proc=0,k=0; proc<size; proc++) { 350 if (!len_s[proc]) continue; 351 /* form outgoing message for i-structure: 352 buf_si[0]: nrows to be sent 353 [1:nrows]: row index (global) 354 [nrows+1:2*nrows+1]: i-structure index 355 */ 356 /*-------------------------------------------*/ 357 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 358 buf_si_i = buf_si + nrows+1; 359 buf_si[0] = nrows; 360 buf_si_i[0] = 0; 361 nrows = 0; 362 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 363 nzi = coi[i+1] - coi[i]; 364 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 365 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 366 nrows++; 367 } 368 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 369 k++; 370 buf_si += len_si[proc]; 371 } 372 for (i=0; i<nrecv; i++) { 373 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 374 } 375 ierr = PetscFree(rwaits);CHKERRQ(ierr); 376 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 377 378 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 379 ierr = PetscFree(len_ri);CHKERRQ(ierr); 380 ierr = PetscFree(swaits);CHKERRQ(ierr); 381 ierr = PetscFree(buf_s);CHKERRQ(ierr); 382 #if defined(PTAP_PROFILE) 383 ierr = PetscTime(&t2);CHKERRQ(ierr); 384 #endif 385 /* (5) compute the local portion of Cmpi */ 386 /* ------------------------------------------ */ 387 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 388 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 389 current_space = free_space; 390 391 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 392 for (k=0; k<nrecv; k++) { 393 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 394 nrows = *buf_ri_k[k]; 395 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 396 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 397 } 398 399 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 400 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 401 for (i=0; i<pn; i++) { 402 /* add C_loc into Cmpi */ 403 nzi = c_loc->i[i+1] - c_loc->i[i]; 404 Jptr = c_loc->j + c_loc->i[i]; 405 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 406 407 /* add received col data into lnk */ 408 for (k=0; k<nrecv; k++) { /* k-th received message */ 409 if (i == *nextrow[k]) { /* i-th row */ 410 nzi = *(nextci[k]+1) - *nextci[k]; 411 Jptr = buf_rj[k] + *nextci[k]; 412 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 413 nextrow[k]++; nextci[k]++; 414 } 415 } 416 nzi = lnk[0]; 417 418 /* copy data into free space, then initialize lnk */ 419 ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 420 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 421 } 422 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 423 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 424 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 425 #if defined(PTAP_PROFILE) 426 ierr = PetscTime(&t3);CHKERRQ(ierr); 427 #endif 428 429 /* (6) create symbolic parallel matrix Cmpi */ 430 /*------------------------------------------*/ 431 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 432 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 433 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 434 ierr = MatSetType(Cmpi,MATMPIAIJ);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 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 454 ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr); 455 ierr = PetscOptionsEnd();CHKERRQ(ierr); 456 457 if (alg == 1) { 458 /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 459 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 460 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 461 } else { 462 Cmpi->ops->ptapnumeric = 0; /* not done yet */ 463 SETERRQ(comm,PETSC_ERR_ARG_SIZ,"Not done yet"); 464 } 465 466 467 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 468 Cmpi->assembled = PETSC_FALSE; 469 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 470 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 471 *C = Cmpi; 472 473 #if defined(PTAP_PROFILE) 474 ierr = PetscTime(&t4);CHKERRQ(ierr); 475 if (rank == 1) { 476 printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0); 477 } 478 #endif 479 PetscFunctionReturn(0); 480 } 481 482 #undef __FUNCT__ 483 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 484 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 485 { 486 PetscErrorCode ierr; 487 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 488 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 489 Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 490 Mat_PtAPMPI *ptap = c->ptap; 491 Mat AP_loc,C_loc,C_oth; 492 PetscInt i,rstart,rend,cm,ncols,row; 493 PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 494 PetscScalar *apa; 495 const PetscInt *cols; 496 const PetscScalar *vals; 497 #if defined(PTAP_PROFILE) 498 PetscMPIInt rank; 499 MPI_Comm comm; 500 PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 501 #endif 502 503 PetscFunctionBegin; 504 ierr = MatZeroEntries(C);CHKERRQ(ierr); 505 506 /* 1) get R = Pd^T,Ro = Po^T */ 507 #if defined(PTAP_PROFILE) 508 ierr = PetscTime(&t0);CHKERRQ(ierr); 509 #endif 510 if (ptap->reuse == MAT_REUSE_MATRIX) { 511 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 512 ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 513 } 514 #if defined(PTAP_PROFILE) 515 ierr = PetscTime(&t1);CHKERRQ(ierr); 516 eR = t1 - t0; 517 #endif 518 519 /* 2) get AP_loc */ 520 AP_loc = ptap->AP_loc; 521 ap = (Mat_SeqAIJ*)AP_loc->data; 522 523 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 524 /*-----------------------------------------------------*/ 525 if (ptap->reuse == MAT_REUSE_MATRIX) { 526 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 527 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 528 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 529 } 530 531 /* 2-2) compute numeric A_loc*P - dominating part */ 532 /* ---------------------------------------------- */ 533 /* get data from symbolic products */ 534 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 535 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 536 apa = ptap->apa; 537 api = ap->i; 538 apj = ap->j; 539 for (i=0; i<am; i++) { 540 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 541 AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 542 apnz = api[i+1] - api[i]; 543 for (j=0; j<apnz; j++) { 544 col = apj[j+api[i]]; 545 ap->a[j+ap->i[i]] = apa[col]; 546 apa[col] = 0.0; 547 } 548 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 549 } 550 #if defined(PTAP_PROFILE) 551 ierr = PetscTime(&t2);CHKERRQ(ierr); 552 eAP = t2 - t1; 553 #endif 554 555 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 556 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 557 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 558 C_loc = ptap->C_loc; 559 C_oth = ptap->C_oth; 560 #if defined(PTAP_PROFILE) 561 ierr = PetscTime(&t3);CHKERRQ(ierr); 562 eCseq = t3 - t2; 563 #endif 564 565 /* add C_loc and Co to to C */ 566 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 567 568 /* C_loc -> C */ 569 cm = C_loc->rmap->N; 570 c_seq = (Mat_SeqAIJ*)C_loc->data; 571 cols = c_seq->j; 572 vals = c_seq->a; 573 for (i=0; i<cm; i++) { 574 ncols = c_seq->i[i+1] - c_seq->i[i]; 575 row = rstart + i; 576 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 577 cols += ncols; vals += ncols; 578 } 579 580 /* Co -> C, off-processor part */ 581 cm = C_oth->rmap->N; 582 c_seq = (Mat_SeqAIJ*)C_oth->data; 583 cols = c_seq->j; 584 vals = c_seq->a; 585 for (i=0; i<cm; i++) { 586 ncols = c_seq->i[i+1] - c_seq->i[i]; 587 row = p->garray[i]; 588 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 589 cols += ncols; vals += ncols; 590 } 591 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 592 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 593 #if defined(PTAP_PROFILE) 594 ierr = PetscTime(&t4);CHKERRQ(ierr); 595 eCmpi = t4 - t3; 596 597 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 598 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 599 if (rank==1) { 600 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); 601 } 602 #endif 603 ptap->reuse = MAT_REUSE_MATRIX; 604 PetscFunctionReturn(0); 605 } 606 607 #undef __FUNCT__ 608 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap" 609 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C) 610 { 611 PetscErrorCode ierr; 612 Mat Cmpi; 613 Mat_PtAPMPI *ptap; 614 PetscFreeSpaceList free_space=NULL,current_space=NULL; 615 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 616 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 617 Mat_SeqAIJ *p_loc,*p_oth; 618 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 619 PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 620 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 621 PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 622 PetscBT lnkbt; 623 MPI_Comm comm; 624 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 625 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 626 PetscInt len,proc,*dnz,*onz,*owners; 627 PetscInt nzi,*pti,*ptj; 628 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 629 MPI_Request *swaits,*rwaits; 630 MPI_Status *sstatus,rstatus; 631 Mat_Merge_SeqsToMPI *merge; 632 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 633 PetscReal afill=1.0,afill_tmp; 634 635 PetscFunctionBegin; 636 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 637 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 638 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 639 640 /* create struct Mat_PtAPMPI and attached it to C later */ 641 ierr = PetscNew(&ptap);CHKERRQ(ierr); 642 ierr = PetscNew(&merge);CHKERRQ(ierr); 643 ptap->merge = merge; 644 ptap->reuse = MAT_INITIAL_MATRIX; 645 646 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 647 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 648 649 /* get P_loc by taking all local rows of P */ 650 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 651 652 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 653 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 654 pi_loc = p_loc->i; pj_loc = p_loc->j; 655 pi_oth = p_oth->i; pj_oth = p_oth->j; 656 657 /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 658 /*--------------------------------------------------------------------------*/ 659 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 660 api[0] = 0; 661 662 /* create and initialize a linked list */ 663 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 664 665 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 666 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 667 current_space = free_space; 668 669 for (i=0; i<am; i++) { 670 /* diagonal portion of A */ 671 nzi = adi[i+1] - adi[i]; 672 aj = ad->j + adi[i]; 673 for (j=0; j<nzi; j++) { 674 row = aj[j]; 675 pnz = pi_loc[row+1] - pi_loc[row]; 676 Jptr = pj_loc + pi_loc[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 of A */ 681 nzi = aoi[i+1] - aoi[i]; 682 aj = ao->j + aoi[i]; 683 for (j=0; j<nzi; j++) { 684 row = aj[j]; 685 pnz = pi_oth[row+1] - pi_oth[row]; 686 Jptr = pj_oth + pi_oth[row]; 687 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 688 } 689 apnz = lnk[0]; 690 api[i+1] = api[i] + apnz; 691 if (ap_rmax < apnz) ap_rmax = apnz; 692 693 /* if free space is not available, double the total space in the list */ 694 if (current_space->local_remaining<apnz) { 695 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 696 nspacedouble++; 697 } 698 699 /* Copy data into free space, then initialize lnk */ 700 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 701 702 current_space->array += apnz; 703 current_space->local_used += apnz; 704 current_space->local_remaining -= apnz; 705 } 706 707 /* Allocate space for apj, initialize apj, and */ 708 /* destroy list of free space and other temporary array(s) */ 709 ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 710 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 711 afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 712 if (afill_tmp > afill) afill = afill_tmp; 713 714 /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 715 /*-----------------------------------------------------------------*/ 716 ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 717 718 /* then, compute symbolic Co = (p->B)^T*AP */ 719 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 720 >= (num of nonzero rows of C_seq) - pn */ 721 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 722 coi[0] = 0; 723 724 /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 725 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am])); 726 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 727 current_space = free_space; 728 729 for (i=0; i<pon; i++) { 730 pnz = poti[i+1] - poti[i]; 731 ptJ = potj + poti[i]; 732 for (j=0; j<pnz; j++) { 733 row = ptJ[j]; /* row of AP == col of Pot */ 734 apnz = api[row+1] - api[row]; 735 Jptr = apj + api[row]; 736 /* add non-zero cols of AP into the sorted linked list lnk */ 737 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 738 } 739 nnz = lnk[0]; 740 741 /* If free space is not available, double the total space in the list */ 742 if (current_space->local_remaining<nnz) { 743 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 744 nspacedouble++; 745 } 746 747 /* Copy data into free space, and zero out denserows */ 748 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 749 750 current_space->array += nnz; 751 current_space->local_used += nnz; 752 current_space->local_remaining -= nnz; 753 754 coi[i+1] = coi[i] + nnz; 755 } 756 757 ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 758 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 759 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 760 if (afill_tmp > afill) afill = afill_tmp; 761 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 762 763 /* (3) send j-array (coj) of Co to other processors */ 764 /*--------------------------------------------------*/ 765 ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 766 len_s = merge->len_s; 767 merge->nsend = 0; 768 769 770 /* determine row ownership */ 771 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 772 merge->rowmap->n = pn; 773 merge->rowmap->bs = 1; 774 775 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 776 owners = merge->rowmap->range; 777 778 /* determine the number of messages to send, their lengths */ 779 ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 780 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 781 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 782 783 proc = 0; 784 for (i=0; i<pon; i++) { 785 while (prmap[i] >= owners[proc+1]) proc++; 786 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 787 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 788 } 789 790 len = 0; /* max length of buf_si[], see (4) */ 791 owners_co[0] = 0; 792 for (proc=0; proc<size; proc++) { 793 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 794 if (len_s[proc]) { 795 merge->nsend++; 796 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 797 len += len_si[proc]; 798 } 799 } 800 801 /* determine the number and length of messages to receive for coi and coj */ 802 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 803 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 804 805 /* post the Irecv and Isend of coj */ 806 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 807 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 808 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 809 for (proc=0, k=0; proc<size; proc++) { 810 if (!len_s[proc]) continue; 811 i = owners_co[proc]; 812 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 813 k++; 814 } 815 816 /* receives and sends of coj are complete */ 817 for (i=0; i<merge->nrecv; i++) { 818 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 819 } 820 ierr = PetscFree(rwaits);CHKERRQ(ierr); 821 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 822 823 /* (4) send and recv coi */ 824 /*-----------------------*/ 825 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 826 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 827 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 828 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 829 for (proc=0,k=0; proc<size; proc++) { 830 if (!len_s[proc]) continue; 831 /* form outgoing message for i-structure: 832 buf_si[0]: nrows to be sent 833 [1:nrows]: row index (global) 834 [nrows+1:2*nrows+1]: i-structure index 835 */ 836 /*-------------------------------------------*/ 837 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 838 buf_si_i = buf_si + nrows+1; 839 buf_si[0] = nrows; 840 buf_si_i[0] = 0; 841 nrows = 0; 842 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 843 nzi = coi[i+1] - coi[i]; 844 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 845 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 846 nrows++; 847 } 848 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 849 k++; 850 buf_si += len_si[proc]; 851 } 852 i = merge->nrecv; 853 while (i--) { 854 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 855 } 856 ierr = PetscFree(rwaits);CHKERRQ(ierr); 857 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 858 859 ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 860 ierr = PetscFree(len_ri);CHKERRQ(ierr); 861 ierr = PetscFree(swaits);CHKERRQ(ierr); 862 ierr = PetscFree(buf_s);CHKERRQ(ierr); 863 864 /* (5) compute the local portion of C (mpi mat) */ 865 /*----------------------------------------------*/ 866 ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 867 868 /* allocate pti array and free space for accumulating nonzero column info */ 869 ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 870 pti[0] = 0; 871 872 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 873 nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am])); 874 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 875 current_space = free_space; 876 877 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 878 for (k=0; k<merge->nrecv; k++) { 879 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 880 nrows = *buf_ri_k[k]; 881 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 882 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 883 } 884 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 885 886 for (i=0; i<pn; i++) { 887 /* add pdt[i,:]*AP into lnk */ 888 pnz = pdti[i+1] - pdti[i]; 889 ptJ = pdtj + pdti[i]; 890 for (j=0; j<pnz; j++) { 891 row = ptJ[j]; /* row of AP == col of Pt */ 892 apnz = api[row+1] - api[row]; 893 Jptr = apj + api[row]; 894 /* add non-zero cols of AP into the sorted linked list lnk */ 895 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 896 } 897 898 /* add received col data into lnk */ 899 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 900 if (i == *nextrow[k]) { /* i-th row */ 901 nzi = *(nextci[k]+1) - *nextci[k]; 902 Jptr = buf_rj[k] + *nextci[k]; 903 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 904 nextrow[k]++; nextci[k]++; 905 } 906 } 907 nnz = lnk[0]; 908 909 /* if free space is not available, make more free space */ 910 if (current_space->local_remaining<nnz) { 911 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 912 nspacedouble++; 913 } 914 /* copy data into free space, then initialize lnk */ 915 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 916 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 917 918 current_space->array += nnz; 919 current_space->local_used += nnz; 920 current_space->local_remaining -= nnz; 921 922 pti[i+1] = pti[i] + nnz; 923 } 924 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 925 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 926 927 ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 928 ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 929 afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 930 if (afill_tmp > afill) afill = afill_tmp; 931 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 932 933 /* (6) create symbolic parallel matrix Cmpi */ 934 /*------------------------------------------*/ 935 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 936 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 937 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 938 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 939 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 940 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 941 942 merge->bi = pti; /* Cseq->i */ 943 merge->bj = ptj; /* Cseq->j */ 944 merge->coi = coi; /* Co->i */ 945 merge->coj = coj; /* Co->j */ 946 merge->buf_ri = buf_ri; 947 merge->buf_rj = buf_rj; 948 merge->owners_co = owners_co; 949 merge->destroy = Cmpi->ops->destroy; 950 951 /* attach the supporting struct to Cmpi for reuse */ 952 c = (Mat_MPIAIJ*)Cmpi->data; 953 c->ptap = ptap; 954 ptap->api = api; 955 ptap->apj = apj; 956 ptap->duplicate = Cmpi->ops->duplicate; 957 ptap->destroy = Cmpi->ops->destroy; 958 959 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 960 Cmpi->assembled = PETSC_FALSE; 961 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 962 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 963 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap; 964 *C = Cmpi; 965 966 /* flag 'scalable' determines which implementations to be used: 967 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 968 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 969 /* set default scalable */ 970 ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */ 971 972 ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 973 if (!ptap->scalable) { /* Do dense axpy */ 974 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 975 } else { 976 ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 977 } 978 979 #if defined(PETSC_USE_INFO) 980 if (pti[pn] != 0) { 981 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 982 ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 983 } else { 984 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 985 } 986 #endif 987 PetscFunctionReturn(0); 988 } 989 990 #undef __FUNCT__ 991 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap" 992 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C) 993 { 994 PetscErrorCode ierr; 995 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 996 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 997 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 998 Mat_SeqAIJ *p_loc,*p_oth; 999 Mat_PtAPMPI *ptap; 1000 Mat_Merge_SeqsToMPI *merge; 1001 PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 1002 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 1003 PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1004 MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1005 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1006 MPI_Comm comm; 1007 PetscMPIInt size,rank,taga,*len_s; 1008 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1009 PetscInt **buf_ri,**buf_rj; 1010 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1011 MPI_Request *s_waits,*r_waits; 1012 MPI_Status *status; 1013 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1014 PetscInt *api,*apj,*coi,*coj; 1015 PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 1016 PetscBool scalable; 1017 #if defined(PTAP_PROFILE) 1018 PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 1019 #endif 1020 1021 PetscFunctionBegin; 1022 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1023 #if defined(PTAP_PROFILE) 1024 ierr = PetscTime(&t0);CHKERRQ(ierr); 1025 #endif 1026 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1027 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1028 1029 ptap = c->ptap; 1030 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"); 1031 merge = ptap->merge; 1032 apa = ptap->apa; 1033 scalable = ptap->scalable; 1034 1035 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 1036 /*-----------------------------------------------------*/ 1037 if (ptap->reuse == MAT_INITIAL_MATRIX) { 1038 /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1039 ptap->reuse = MAT_REUSE_MATRIX; 1040 } else { /* update numerical values of P_oth and P_loc */ 1041 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1042 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1043 } 1044 #if defined(PTAP_PROFILE) 1045 ierr = PetscTime(&t1);CHKERRQ(ierr); 1046 eP = t1-t0; 1047 #endif 1048 /* 1049 printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 1050 a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 1051 ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 1052 ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 1053 */ 1054 1055 /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1056 /*--------------------------------------------------------------*/ 1057 /* get data from symbolic products */ 1058 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1059 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1060 pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 1061 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 1062 1063 coi = merge->coi; coj = merge->coj; 1064 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1065 1066 bi = merge->bi; bj = merge->bj; 1067 owners = merge->rowmap->range; 1068 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 1069 1070 api = ptap->api; apj = ptap->apj; 1071 1072 if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1073 ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 1074 #if defined(PTAP_PROFILE) 1075 ierr = PetscTime(&t1);CHKERRQ(ierr); 1076 #endif 1077 for (i=0; i<am; i++) { 1078 #if defined(PTAP_PROFILE) 1079 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1080 #endif 1081 /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1082 /*------------------------------------------------------------*/ 1083 apJ = apj + api[i]; 1084 1085 /* diagonal portion of A */ 1086 anz = adi[i+1] - adi[i]; 1087 adj = ad->j + adi[i]; 1088 ada = ad->a + adi[i]; 1089 for (j=0; j<anz; j++) { 1090 row = adj[j]; 1091 pnz = pi_loc[row+1] - pi_loc[row]; 1092 pj = pj_loc + pi_loc[row]; 1093 pa = pa_loc + pi_loc[row]; 1094 1095 /* perform dense axpy */ 1096 valtmp = ada[j]; 1097 for (k=0; k<pnz; k++) { 1098 apa[pj[k]] += valtmp*pa[k]; 1099 } 1100 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1101 } 1102 1103 /* off-diagonal portion of A */ 1104 anz = aoi[i+1] - aoi[i]; 1105 aoj = ao->j + aoi[i]; 1106 aoa = ao->a + aoi[i]; 1107 for (j=0; j<anz; j++) { 1108 row = aoj[j]; 1109 pnz = pi_oth[row+1] - pi_oth[row]; 1110 pj = pj_oth + pi_oth[row]; 1111 pa = pa_oth + pi_oth[row]; 1112 1113 /* perform dense axpy */ 1114 valtmp = aoa[j]; 1115 for (k=0; k<pnz; k++) { 1116 apa[pj[k]] += valtmp*pa[k]; 1117 } 1118 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1119 } 1120 #if defined(PTAP_PROFILE) 1121 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1122 et2_AP += t2_1 - t2_0; 1123 #endif 1124 1125 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1126 /*--------------------------------------------------------------*/ 1127 apnz = api[i+1] - api[i]; 1128 /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1129 pnz = po->i[i+1] - po->i[i]; 1130 poJ = po->j + po->i[i]; 1131 pA = po->a + po->i[i]; 1132 for (j=0; j<pnz; j++) { 1133 row = poJ[j]; 1134 cnz = coi[row+1] - coi[row]; 1135 cj = coj + coi[row]; 1136 ca = coa + coi[row]; 1137 /* perform dense axpy */ 1138 valtmp = pA[j]; 1139 for (k=0; k<cnz; k++) { 1140 ca[k] += valtmp*apa[cj[k]]; 1141 } 1142 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1143 } 1144 /* put the value into Cd (diagonal part) */ 1145 pnz = pd->i[i+1] - pd->i[i]; 1146 pdJ = pd->j + pd->i[i]; 1147 pA = pd->a + pd->i[i]; 1148 for (j=0; j<pnz; j++) { 1149 row = pdJ[j]; 1150 cnz = bi[row+1] - bi[row]; 1151 cj = bj + bi[row]; 1152 ca = ba + bi[row]; 1153 /* perform dense axpy */ 1154 valtmp = pA[j]; 1155 for (k=0; k<cnz; k++) { 1156 ca[k] += valtmp*apa[cj[k]]; 1157 } 1158 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1159 } 1160 1161 /* zero the current row of A*P */ 1162 for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1163 #if defined(PTAP_PROFILE) 1164 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 1165 ePtAP += t2_2 - t2_1; 1166 #endif 1167 } 1168 } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1169 ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 1170 /*-----------------------------------------------------------------------------------------*/ 1171 pA=pa_loc; 1172 for (i=0; i<am; i++) { 1173 #if defined(PTAP_PROFILE) 1174 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1175 #endif 1176 /* form i-th sparse row of A*P */ 1177 apnz = api[i+1] - api[i]; 1178 apJ = apj + api[i]; 1179 /* diagonal portion of A */ 1180 anz = adi[i+1] - adi[i]; 1181 adj = ad->j + adi[i]; 1182 ada = ad->a + adi[i]; 1183 for (j=0; j<anz; j++) { 1184 row = adj[j]; 1185 pnz = pi_loc[row+1] - pi_loc[row]; 1186 pj = pj_loc + pi_loc[row]; 1187 pa = pa_loc + pi_loc[row]; 1188 valtmp = ada[j]; 1189 nextp = 0; 1190 for (k=0; nextp<pnz; k++) { 1191 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1192 apa[k] += valtmp*pa[nextp++]; 1193 } 1194 } 1195 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1196 } 1197 /* off-diagonal portion of A */ 1198 anz = aoi[i+1] - aoi[i]; 1199 aoj = ao->j + aoi[i]; 1200 aoa = ao->a + aoi[i]; 1201 for (j=0; j<anz; j++) { 1202 row = aoj[j]; 1203 pnz = pi_oth[row+1] - pi_oth[row]; 1204 pj = pj_oth + pi_oth[row]; 1205 pa = pa_oth + pi_oth[row]; 1206 valtmp = aoa[j]; 1207 nextp = 0; 1208 for (k=0; nextp<pnz; k++) { 1209 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1210 apa[k] += valtmp*pa[nextp++]; 1211 } 1212 } 1213 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1214 } 1215 #if defined(PTAP_PROFILE) 1216 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1217 et2_AP += t2_1 - t2_0; 1218 #endif 1219 1220 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1221 /*--------------------------------------------------------------*/ 1222 pnz = pi_loc[i+1] - pi_loc[i]; 1223 pJ = pj_loc + pi_loc[i]; 1224 for (j=0; j<pnz; j++) { 1225 nextap = 0; 1226 row = pJ[j]; /* global index */ 1227 if (row < pcstart || row >=pcend) { /* put the value into Co */ 1228 row = *poJ; 1229 cj = coj + coi[row]; 1230 ca = coa + coi[row]; poJ++; 1231 } else { /* put the value into Cd */ 1232 row = *pdJ; 1233 cj = bj + bi[row]; 1234 ca = ba + bi[row]; pdJ++; 1235 } 1236 valtmp = pA[j]; 1237 for (k=0; nextap<apnz; k++) { 1238 if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 1239 } 1240 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1241 } 1242 pA += pnz; 1243 /* zero the current row info for A*P */ 1244 ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1245 #if defined(PTAP_PROFILE) 1246 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 1247 ePtAP += t2_2 - t2_1; 1248 #endif 1249 } 1250 } 1251 #if defined(PTAP_PROFILE) 1252 ierr = PetscTime(&t2);CHKERRQ(ierr); 1253 #endif 1254 1255 /* 3) send and recv matrix values coa */ 1256 /*------------------------------------*/ 1257 buf_ri = merge->buf_ri; 1258 buf_rj = merge->buf_rj; 1259 len_s = merge->len_s; 1260 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1261 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1262 1263 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1264 for (proc=0,k=0; proc<size; proc++) { 1265 if (!len_s[proc]) continue; 1266 i = merge->owners_co[proc]; 1267 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1268 k++; 1269 } 1270 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1271 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1272 1273 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1274 ierr = PetscFree(r_waits);CHKERRQ(ierr); 1275 ierr = PetscFree(coa);CHKERRQ(ierr); 1276 #if defined(PTAP_PROFILE) 1277 ierr = PetscTime(&t3);CHKERRQ(ierr); 1278 #endif 1279 1280 /* 4) insert local Cseq and received values into Cmpi */ 1281 /*------------------------------------------------------*/ 1282 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1283 for (k=0; k<merge->nrecv; k++) { 1284 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1285 nrows = *(buf_ri_k[k]); 1286 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1287 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1288 } 1289 1290 for (i=0; i<cm; i++) { 1291 row = owners[rank] + i; /* global row index of C_seq */ 1292 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1293 ba_i = ba + bi[i]; 1294 bnz = bi[i+1] - bi[i]; 1295 /* add received vals into ba */ 1296 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1297 /* i-th row */ 1298 if (i == *nextrow[k]) { 1299 cnz = *(nextci[k]+1) - *nextci[k]; 1300 cj = buf_rj[k] + *(nextci[k]); 1301 ca = abuf_r[k] + *(nextci[k]); 1302 nextcj = 0; 1303 for (j=0; nextcj<cnz; j++) { 1304 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1305 ba_i[j] += ca[nextcj++]; 1306 } 1307 } 1308 nextrow[k]++; nextci[k]++; 1309 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1310 } 1311 } 1312 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1313 } 1314 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1315 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1316 1317 ierr = PetscFree(ba);CHKERRQ(ierr); 1318 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1319 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1320 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1321 #if defined(PTAP_PROFILE) 1322 ierr = PetscTime(&t4);CHKERRQ(ierr); 1323 if (rank==1) { 1324 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); 1325 } 1326 #endif 1327 PetscFunctionReturn(0); 1328 } 1329