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