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