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