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