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