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