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