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 #include <petsc/private/hashmapiv.h> 13 #include <petsc/private/hashseti.h> 14 #include <petscsf.h> 15 16 17 PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer) 18 { 19 PetscErrorCode ierr; 20 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 21 Mat_APMPI *ptap=a->ap; 22 PetscBool iascii; 23 PetscViewerFormat format; 24 25 PetscFunctionBegin; 26 if (!ptap) { 27 /* hack: MatDuplicate() sets oldmat->ops->view to newmat which is a base mat class with null ptpa! */ 28 A->ops->view = MatView_MPIAIJ; 29 ierr = (A->ops->view)(A,viewer);CHKERRQ(ierr); 30 PetscFunctionReturn(0); 31 } 32 33 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 34 if (iascii) { 35 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 36 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 37 if (ptap->algType == 0) { 38 ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr); 39 } else if (ptap->algType == 1) { 40 ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr); 41 } else if (ptap->algType == 2) { 42 ierr = PetscViewerASCIIPrintf(viewer,"using allatonce MatPtAP() implementation\n");CHKERRQ(ierr); 43 } else if (ptap->algType == 3) { 44 ierr = PetscViewerASCIIPrintf(viewer,"using merged allatonce MatPtAP() implementation\n");CHKERRQ(ierr); 45 } 46 } 47 } 48 ierr = (ptap->view)(A,viewer);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 52 PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_AP(Mat A) 53 { 54 PetscErrorCode ierr; 55 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 56 Mat_APMPI *ptap=a->ap; 57 Mat_Merge_SeqsToMPI *merge; 58 59 PetscFunctionBegin; 60 if (!ptap) PetscFunctionReturn(0); 61 62 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 63 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 64 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 65 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 66 ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 67 ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 68 ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 69 if (ptap->AP_loc) { /* used by alg_rap */ 70 Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 71 ierr = PetscFree(ap->i);CHKERRQ(ierr); 72 ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 73 ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 74 } else { /* used by alg_ptap */ 75 ierr = PetscFree(ptap->api);CHKERRQ(ierr); 76 ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 77 } 78 ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 79 ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 80 if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 81 82 ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr); 83 84 merge=ptap->merge; 85 if (merge) { /* used by alg_ptap */ 86 ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 87 ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 88 ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 89 ierr = PetscFree(merge->bi);CHKERRQ(ierr); 90 ierr = PetscFree(merge->bj);CHKERRQ(ierr); 91 ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 92 ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 93 ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 94 ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 95 ierr = PetscFree(merge->coi);CHKERRQ(ierr); 96 ierr = PetscFree(merge->coj);CHKERRQ(ierr); 97 ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 98 ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 99 ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 100 } 101 ierr = ISLocalToGlobalMappingDestroy(&ptap->ltog);CHKERRQ(ierr); 102 103 ierr = PetscSFDestroy(&ptap->sf);CHKERRQ(ierr); 104 ierr = PetscFree(ptap->c_othi);CHKERRQ(ierr); 105 ierr = PetscFree(ptap->c_rmti);CHKERRQ(ierr); 106 PetscFunctionReturn(0); 107 } 108 109 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 110 { 111 PetscErrorCode ierr; 112 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 113 Mat_APMPI *ptap=a->ap; 114 115 PetscFunctionBegin; 116 ierr = (*A->ops->freeintermediatedatastructures)(A);CHKERRQ(ierr); 117 ierr = (*ptap->destroy)(A);CHKERRQ(ierr); /* MatDestroy_MPIAIJ(A) */ 118 ierr = PetscFree(ptap);CHKERRQ(ierr); 119 PetscFunctionReturn(0); 120 } 121 122 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 123 { 124 PetscErrorCode ierr; 125 PetscBool flg; 126 MPI_Comm comm; 127 #if !defined(PETSC_HAVE_HYPRE) 128 const char *algTypes[4] = {"scalable","nonscalable","allatonce","allatonce_merged"}; 129 PetscInt nalg=4; 130 #else 131 const char *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","hypre"}; 132 PetscInt nalg=5; 133 #endif 134 PetscInt pN=P->cmap->N,alg=1; /* set default algorithm */ 135 136 PetscFunctionBegin; 137 /* check if matrix local sizes are compatible */ 138 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 139 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); 140 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); 141 142 if (scall == MAT_INITIAL_MATRIX) { 143 /* pick an algorithm */ 144 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 145 ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 146 ierr = PetscOptionsEnd();CHKERRQ(ierr); 147 148 if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */ 149 MatInfo Ainfo,Pinfo; 150 PetscInt nz_local; 151 PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 152 153 ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 154 ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 155 nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); 156 157 if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 158 ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 159 160 if (alg_scalable) { 161 alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 162 } 163 } 164 165 switch (alg) { 166 case 1: 167 /* do R=P^T locally, then C=R*A*P -- nonscalable */ 168 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 169 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 170 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 171 break; 172 case 2: 173 /* compute C=P^T*A*P allatonce */ 174 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 175 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(A,P,fill,C);CHKERRQ(ierr); 176 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 177 break; 178 case 3: 179 /* compute C=P^T*A*P allatonce */ 180 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 181 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(A,P,fill,C);CHKERRQ(ierr); 182 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 183 break; 184 #if defined(PETSC_HAVE_HYPRE) 185 case 4: 186 /* Use boomerAMGBuildCoarseOperator */ 187 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 188 ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 189 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 190 break; 191 #endif 192 default: 193 /* do R=P^T locally, then C=R*A*P */ 194 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 195 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr); 196 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 197 break; 198 } 199 200 if (alg == 0 || alg == 1 || alg == 2 || alg == 3) { 201 Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 202 Mat_APMPI *ap = c->ap; 203 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 204 ap->freestruct = PETSC_FALSE; 205 ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 206 ierr = PetscOptionsEnd();CHKERRQ(ierr); 207 } 208 } 209 210 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 211 ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 212 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 213 PetscFunctionReturn(0); 214 } 215 216 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C) 217 { 218 PetscErrorCode ierr; 219 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 220 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 221 Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 222 Mat_APMPI *ptap = c->ap; 223 Mat AP_loc,C_loc,C_oth; 224 PetscInt i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz,nout; 225 PetscScalar *apa; 226 const PetscInt *cols; 227 const PetscScalar *vals; 228 229 PetscFunctionBegin; 230 if (!ptap->AP_loc) { 231 MPI_Comm comm; 232 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 233 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 234 } 235 236 ierr = MatZeroEntries(C);CHKERRQ(ierr); 237 238 /* 1) get R = Pd^T,Ro = Po^T */ 239 if (ptap->reuse == MAT_REUSE_MATRIX) { 240 ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 241 ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 242 } 243 244 /* 2) get AP_loc */ 245 AP_loc = ptap->AP_loc; 246 ap = (Mat_SeqAIJ*)AP_loc->data; 247 248 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 249 /*-----------------------------------------------------*/ 250 if (ptap->reuse == MAT_REUSE_MATRIX) { 251 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 252 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 253 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 254 } 255 256 /* 2-2) compute numeric A_loc*P - dominating part */ 257 /* ---------------------------------------------- */ 258 /* get data from symbolic products */ 259 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 260 if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 261 262 api = ap->i; 263 apj = ap->j; 264 ierr = ISLocalToGlobalMappingApply(ptap->ltog,api[AP_loc->rmap->n],apj,apj);CHKERRQ(ierr); 265 for (i=0; i<am; i++) { 266 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 267 apnz = api[i+1] - api[i]; 268 apa = ap->a + api[i]; 269 ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr); 270 AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa); 271 } 272 ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,api[AP_loc->rmap->n],apj,&nout,apj);CHKERRQ(ierr); 273 if (api[AP_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",api[AP_loc->rmap->n],nout); 274 275 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 276 /* Always use scalable version since we are in the MPI scalable version */ 277 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 278 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 279 280 C_loc = ptap->C_loc; 281 C_oth = ptap->C_oth; 282 283 /* add C_loc and Co to to C */ 284 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 285 286 /* C_loc -> C */ 287 cm = C_loc->rmap->N; 288 c_seq = (Mat_SeqAIJ*)C_loc->data; 289 cols = c_seq->j; 290 vals = c_seq->a; 291 ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_loc->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr); 292 293 /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 294 /* when there are no off-processor parts. */ 295 /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 296 /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 297 /* a table, and other, more complex stuff has to be done. */ 298 if (C->assembled) { 299 C->was_assembled = PETSC_TRUE; 300 C->assembled = PETSC_FALSE; 301 } 302 if (C->was_assembled) { 303 for (i=0; i<cm; i++) { 304 ncols = c_seq->i[i+1] - c_seq->i[i]; 305 row = rstart + i; 306 ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 307 cols += ncols; vals += ncols; 308 } 309 } else { 310 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 311 } 312 ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_loc->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr); 313 if (c_seq->i[C_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_seq->i[C_loc->rmap->n],nout); 314 315 /* Co -> C, off-processor part */ 316 cm = C_oth->rmap->N; 317 c_seq = (Mat_SeqAIJ*)C_oth->data; 318 cols = c_seq->j; 319 vals = c_seq->a; 320 ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_oth->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr); 321 for (i=0; i<cm; i++) { 322 ncols = c_seq->i[i+1] - c_seq->i[i]; 323 row = p->garray[i]; 324 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 325 cols += ncols; vals += ncols; 326 } 327 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 328 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 329 330 ptap->reuse = MAT_REUSE_MATRIX; 331 332 ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_oth->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr); 333 if (c_seq->i[C_oth->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_seq->i[C_loc->rmap->n],nout); 334 335 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 336 if (ptap->freestruct) { 337 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 338 } 339 PetscFunctionReturn(0); 340 } 341 342 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C) 343 { 344 PetscErrorCode ierr; 345 Mat_APMPI *ptap; 346 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 347 MPI_Comm comm; 348 PetscMPIInt size,rank; 349 Mat Cmpi,P_loc,P_oth; 350 PetscFreeSpaceList free_space=NULL,current_space=NULL; 351 PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 352 PetscInt *lnk,i,k,pnz,row,nsend; 353 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 354 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 355 PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 356 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 357 MPI_Request *swaits,*rwaits; 358 MPI_Status *sstatus,rstatus; 359 PetscLayout rowmap; 360 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 361 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 362 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi,nout; 363 Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 364 PetscScalar *apv; 365 PetscTable ta; 366 MatType mtype; 367 const char *prefix; 368 #if defined(PETSC_USE_INFO) 369 PetscReal apfill; 370 #endif 371 372 PetscFunctionBegin; 373 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 374 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 375 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 376 377 if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 378 379 /* create symbolic parallel matrix Cmpi */ 380 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 381 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 382 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 383 384 /* create struct Mat_APMPI and attached it to C later */ 385 ierr = PetscNew(&ptap);CHKERRQ(ierr); 386 ptap->reuse = MAT_INITIAL_MATRIX; 387 ptap->algType = 0; 388 389 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 390 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr); 391 /* get P_loc by taking all local rows of P */ 392 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr); 393 394 ptap->P_loc = P_loc; 395 ptap->P_oth = P_oth; 396 397 /* (0) compute Rd = Pd^T, Ro = Po^T */ 398 /* --------------------------------- */ 399 ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 400 ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 401 402 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 403 /* ----------------------------------------------------------------- */ 404 p_loc = (Mat_SeqAIJ*)P_loc->data; 405 if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data; 406 407 /* create and initialize a linked list */ 408 ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 409 MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta); 410 MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta); 411 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 412 413 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 414 415 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 416 if (ao) { 417 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 418 } else { 419 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 420 } 421 current_space = free_space; 422 nspacedouble = 0; 423 424 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 425 api[0] = 0; 426 for (i=0; i<am; i++) { 427 /* diagonal portion: Ad[i,:]*P */ 428 ai = ad->i; pi = p_loc->i; 429 nzi = ai[i+1] - ai[i]; 430 aj = ad->j + ai[i]; 431 for (j=0; j<nzi; j++) { 432 row = aj[j]; 433 pnz = pi[row+1] - pi[row]; 434 Jptr = p_loc->j + pi[row]; 435 /* add non-zero cols of P into the sorted linked list lnk */ 436 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 437 } 438 /* off-diagonal portion: Ao[i,:]*P */ 439 if (ao) { 440 ai = ao->i; pi = p_oth->i; 441 nzi = ai[i+1] - ai[i]; 442 aj = ao->j + ai[i]; 443 for (j=0; j<nzi; j++) { 444 row = aj[j]; 445 pnz = pi[row+1] - pi[row]; 446 Jptr = p_oth->j + pi[row]; 447 ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 448 } 449 } 450 apnz = lnk[0]; 451 api[i+1] = api[i] + apnz; 452 453 /* if free space is not available, double the total space in the list */ 454 if (current_space->local_remaining<apnz) { 455 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 456 nspacedouble++; 457 } 458 459 /* Copy data into free space, then initialize lnk */ 460 ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 461 462 current_space->array += apnz; 463 current_space->local_used += apnz; 464 current_space->local_remaining -= apnz; 465 } 466 /* Allocate space for apj and apv, initialize apj, and */ 467 /* destroy list of free space and other temporary array(s) */ 468 ierr = PetscCalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 469 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 470 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 471 472 /* Create AP_loc for reuse */ 473 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 474 ierr = MatSeqAIJCompactOutExtraColumns_SeqAIJ(ptap->AP_loc, &ptap->ltog);CHKERRQ(ierr); 475 476 #if defined(PETSC_USE_INFO) 477 if (ao) { 478 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 479 } else { 480 apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 481 } 482 ptap->AP_loc->info.mallocs = nspacedouble; 483 ptap->AP_loc->info.fill_ratio_given = fill; 484 ptap->AP_loc->info.fill_ratio_needed = apfill; 485 486 if (api[am]) { 487 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); 488 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 489 } else { 490 ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 491 } 492 #endif 493 494 /* (2-1) compute symbolic Co = Ro*AP_loc */ 495 /* ------------------------------------ */ 496 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 497 ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 498 ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 499 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 500 501 /* (3) send coj of C_oth to other processors */ 502 /* ------------------------------------------ */ 503 /* determine row ownership */ 504 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 505 rowmap->n = pn; 506 rowmap->bs = 1; 507 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 508 owners = rowmap->range; 509 510 /* determine the number of messages to send, their lengths */ 511 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 512 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 513 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 514 515 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 516 coi = c_oth->i; coj = c_oth->j; 517 con = ptap->C_oth->rmap->n; 518 proc = 0; 519 ierr = ISLocalToGlobalMappingApply(ptap->ltog,coi[con],coj,coj);CHKERRQ(ierr); 520 for (i=0; i<con; i++) { 521 while (prmap[i] >= owners[proc+1]) proc++; 522 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 523 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 524 } 525 526 len = 0; /* max length of buf_si[], see (4) */ 527 owners_co[0] = 0; 528 nsend = 0; 529 for (proc=0; proc<size; proc++) { 530 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 531 if (len_s[proc]) { 532 nsend++; 533 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 534 len += len_si[proc]; 535 } 536 } 537 538 /* determine the number and length of messages to receive for coi and coj */ 539 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 540 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 541 542 /* post the Irecv and Isend of coj */ 543 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 544 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 545 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 546 for (proc=0, k=0; proc<size; proc++) { 547 if (!len_s[proc]) continue; 548 i = owners_co[proc]; 549 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 550 k++; 551 } 552 553 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 554 /* ---------------------------------------- */ 555 ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 556 ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 557 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 558 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 559 ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,c_loc->j);CHKERRQ(ierr); 560 561 /* receives coj are complete */ 562 for (i=0; i<nrecv; i++) { 563 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 564 } 565 ierr = PetscFree(rwaits);CHKERRQ(ierr); 566 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 567 568 /* add received column indices into ta to update Crmax */ 569 for (k=0; k<nrecv; k++) {/* k-th received message */ 570 Jptr = buf_rj[k]; 571 for (j=0; j<len_r[k]; j++) { 572 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 573 } 574 } 575 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 576 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 577 578 /* (4) send and recv coi */ 579 /*-----------------------*/ 580 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 581 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 582 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 583 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 584 for (proc=0,k=0; proc<size; proc++) { 585 if (!len_s[proc]) continue; 586 /* form outgoing message for i-structure: 587 buf_si[0]: nrows to be sent 588 [1:nrows]: row index (global) 589 [nrows+1:2*nrows+1]: i-structure index 590 */ 591 /*-------------------------------------------*/ 592 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 593 buf_si_i = buf_si + nrows+1; 594 buf_si[0] = nrows; 595 buf_si_i[0] = 0; 596 nrows = 0; 597 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 598 nzi = coi[i+1] - coi[i]; 599 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 600 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 601 nrows++; 602 } 603 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 604 k++; 605 buf_si += len_si[proc]; 606 } 607 for (i=0; i<nrecv; i++) { 608 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 609 } 610 ierr = PetscFree(rwaits);CHKERRQ(ierr); 611 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 612 613 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 614 ierr = PetscFree(len_ri);CHKERRQ(ierr); 615 ierr = PetscFree(swaits);CHKERRQ(ierr); 616 ierr = PetscFree(buf_s);CHKERRQ(ierr); 617 618 /* (5) compute the local portion of Cmpi */ 619 /* ------------------------------------------ */ 620 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 621 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 622 current_space = free_space; 623 624 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 625 for (k=0; k<nrecv; k++) { 626 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 627 nrows = *buf_ri_k[k]; 628 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 629 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 630 } 631 632 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 633 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 634 for (i=0; i<pn; i++) { 635 /* add C_loc into Cmpi */ 636 nzi = c_loc->i[i+1] - c_loc->i[i]; 637 Jptr = c_loc->j + c_loc->i[i]; 638 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 639 640 /* add received col data into lnk */ 641 for (k=0; k<nrecv; k++) { /* k-th received message */ 642 if (i == *nextrow[k]) { /* i-th row */ 643 nzi = *(nextci[k]+1) - *nextci[k]; 644 Jptr = buf_rj[k] + *nextci[k]; 645 ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 646 nextrow[k]++; nextci[k]++; 647 } 648 } 649 nzi = lnk[0]; 650 651 /* copy data into free space, then initialize lnk */ 652 ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr); 653 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 654 } 655 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 656 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 657 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 658 659 /* local sizes and preallocation */ 660 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 661 if (P->cmap->bs > 0) { 662 ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr); 663 ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr); 664 } 665 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 666 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 667 668 /* members in merge */ 669 ierr = PetscFree(id_r);CHKERRQ(ierr); 670 ierr = PetscFree(len_r);CHKERRQ(ierr); 671 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 672 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 673 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 674 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 675 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 676 677 /* attach the supporting struct to Cmpi for reuse */ 678 c = (Mat_MPIAIJ*)Cmpi->data; 679 c->ap = ptap; 680 ptap->duplicate = Cmpi->ops->duplicate; 681 ptap->destroy = Cmpi->ops->destroy; 682 ptap->view = Cmpi->ops->view; 683 684 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 685 Cmpi->assembled = PETSC_FALSE; 686 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 687 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 688 Cmpi->ops->view = MatView_MPIAIJ_PtAP; 689 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 690 *C = Cmpi; 691 692 nout = 0; 693 ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_oth->i[ptap->C_oth->rmap->n],c_oth->j,&nout,c_oth->j);CHKERRQ(ierr); 694 if (c_oth->i[ptap->C_oth->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_oth->i[ptap->C_oth->rmap->n],nout); 695 ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,&nout,c_loc->j);CHKERRQ(ierr); 696 if (c_loc->i[ptap->C_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_loc->i[ptap->C_loc->rmap->n],nout); 697 698 PetscFunctionReturn(0); 699 } 700 701 PETSC_STATIC_INLINE PetscErrorCode MatPtAPSymbolicComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,PetscInt i,PetscHSetI dht,PetscHSetI oht) 702 { 703 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 704 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_oth=(Mat_SeqAIJ*)P_oth->data,*pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data; 705 PetscInt *ai,nzi,j,*aj,row,col,*pi,*pj,pnz,nzpi,*p_othcols,k; 706 PetscInt pcstart,pcend,column; 707 PetscErrorCode ierr; 708 709 PetscFunctionBegin; 710 711 pcstart = P->cmap->rstart; 712 pcend = P->cmap->rend; 713 /* diagonal portion: Ad[i,:]*P */ 714 ai = ad->i; 715 nzi = ai[i+1] - ai[i]; 716 aj = ad->j + ai[i]; 717 for (j=0; j<nzi; j++) { 718 row = aj[j]; 719 nzpi = pd->i[row+1] - pd->i[row]; 720 pj = pd->j + pd->i[row]; 721 for (k=0; k<nzpi; k++) { 722 ierr = PetscHSetIAdd(dht,pj[k]+pcstart);CHKERRQ(ierr); 723 } 724 } 725 for (j=0; j<nzi; j++) { 726 row = aj[j]; 727 nzpi = po->i[row+1] - po->i[row]; 728 pj = po->j + po->i[row]; 729 for (k=0; k<nzpi; k++) { 730 ierr = PetscHSetIAdd(oht,p->garray[pj[k]]);CHKERRQ(ierr); 731 } 732 } 733 734 /* off diagonal part: Ao[i, :]*P_oth */ 735 if (ao) { 736 ai = ao->i; 737 pi = p_oth->i; 738 nzi = ai[i+1] - ai[i]; 739 aj = ao->j + ai[i]; 740 for (j=0; j<nzi; j++) { 741 row = aj[j]; 742 pnz = pi[row+1] - pi[row]; 743 p_othcols = p_oth->j + pi[row]; 744 for (col=0; col<pnz; col++) { 745 column = p_othcols[col]; 746 if (column>=pcstart && column<pcend) { 747 ierr = PetscHSetIAdd(dht,column);CHKERRQ(ierr); 748 } else { 749 ierr = PetscHSetIAdd(oht,column);CHKERRQ(ierr); 750 } 751 } 752 } 753 } /* end if (ao) */ 754 PetscFunctionReturn(0); 755 } 756 757 PETSC_STATIC_INLINE PetscErrorCode MatPtAPNumericComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,PetscInt i,PetscHMapIV hmap) 758 { 759 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 760 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_oth=(Mat_SeqAIJ*)P_oth->data,*pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data; 761 PetscInt *ai,nzi,j,*aj,row,col,*pi,pnz,*p_othcols,pcstart,*pj,k,nzpi; 762 PetscScalar ra,*aa,*pa; 763 PetscLogDouble flops=0.0; 764 PetscErrorCode ierr; 765 766 PetscFunctionBegin; 767 pcstart = P->cmap->rstart; 768 /* diagonal portion: Ad[i,:]*P */ 769 ai = ad->i; 770 nzi = ai[i+1] - ai[i]; 771 aj = ad->j + ai[i]; 772 aa = ad->a + ai[i]; 773 774 for (j=0; j<nzi; j++) { 775 ra = aa[j]; 776 row = aj[j]; 777 nzpi = pd->i[row+1] - pd->i[row]; 778 pj = pd->j + pd->i[row]; 779 pa = pd->a + pd->i[row]; 780 for (k=0; k<nzpi; k++) { 781 ierr = PetscHMapIVAddValue(hmap,pj[k]+pcstart,ra*pa[k]);CHKERRQ(ierr); 782 } 783 flops += nzpi; 784 } 785 for (j=0; j<nzi; j++) { 786 ra = aa[j]; 787 row = aj[j]; 788 nzpi = po->i[row+1] - po->i[row]; 789 pj = po->j + po->i[row]; 790 pa = po->a + po->i[row]; 791 for (k=0; k<nzpi; k++) { 792 ierr = PetscHMapIVAddValue(hmap,p->garray[pj[k]],ra*pa[k]);CHKERRQ(ierr); 793 } 794 flops += nzpi; 795 } 796 797 798 /* off diagonal part: Ao[i, :]*P_oth */ 799 if (ao) { 800 ai = ao->i; 801 pi = p_oth->i; 802 nzi = ai[i+1] - ai[i]; 803 aj = ao->j + ai[i]; 804 aa = ao->a + ai[i]; 805 for (j=0; j<nzi; j++) { 806 row = aj[j]; 807 ra = aa[j]; 808 pnz = pi[row+1] - pi[row]; 809 p_othcols = p_oth->j + pi[row]; 810 pa = p_oth->a + pi[row]; 811 for (col=0; col<pnz; col++) { 812 ierr = PetscHMapIVAddValue(hmap,p_othcols[col],ra*pa[col]);CHKERRQ(ierr); 813 } 814 flops += pnz; 815 } 816 } /* end if (ao) */ 817 818 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 819 PetscFunctionReturn(0); 820 } 821 822 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,Mat C) 823 { 824 PetscErrorCode ierr; 825 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 826 Mat_SeqAIJ *cd,*co,*po,*pd; 827 Mat_APMPI *ptap = c->ap; 828 PetscHMapIV hmap; 829 PetscInt i,j,jj,kk,nzi,*c_rmtj,voff,*c_othj,pn,pon,pcstart,pcend,ccstart,ccend,row,am,*poj,*pdj,*apindices,cmaxr,*c_rmtc,*c_rmtjj,*dcc,*occ,loc; 830 PetscScalar *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa; 831 PetscLogDouble flops=0.0; 832 MPI_Comm comm; 833 834 PetscFunctionBegin; 835 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 836 if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 837 838 ierr = MatZeroEntries(C);CHKERRQ(ierr); 839 840 /* Get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 841 /*-----------------------------------------------------*/ 842 if (ptap->reuse == MAT_REUSE_MATRIX) { 843 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 844 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 845 } 846 847 po = (Mat_SeqAIJ*) p->B->data; 848 pd = (Mat_SeqAIJ*) p->A->data; 849 850 ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 851 852 ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr); 853 ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 854 cmaxr = 0; 855 for (i=0; i<pon; i++) { 856 cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]); 857 } 858 ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr); 859 ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 860 ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 861 for (i=0; i<am && pon; i++) { 862 ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 863 nzi = po->i[i+1] - po->i[i]; 864 if (!nzi) continue; 865 ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr); 866 voff = 0; 867 ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 868 if (!voff) continue; 869 /*ierr = PetscMemzero(c_rmtc,sizeof(PetscInt)*pon);CHKERRQ(ierr);*/ 870 871 /* Form C(ii, :) */ 872 poj = po->j + po->i[i]; 873 poa = po->a + po->i[i]; 874 for (j=0; j<nzi; j++) { 875 c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]]; 876 c_rmtaa = c_rmta + ptap->c_rmti[poj[j]]; 877 for (jj=0; jj<voff; jj++) { 878 apvaluestmp[jj] = apvalues[jj]*poa[j]; 879 /*If the row is empty */ 880 if (!c_rmtc[poj[j]]) { 881 c_rmtjj[jj] = apindices[jj]; 882 c_rmtaa[jj] = apvaluestmp[jj]; 883 c_rmtc[poj[j]]++; 884 } else { 885 ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr); 886 if (loc>=0){ /* hit */ 887 c_rmtaa[loc] += apvaluestmp[jj]; 888 } else { /* new element */ 889 loc = -(loc+1); 890 /* Move data backward */ 891 for (kk=c_rmtc[poj[j]]; kk>loc; kk--) { 892 c_rmtjj[kk] = c_rmtjj[kk-1]; 893 c_rmtaa[kk] = c_rmtaa[kk-1]; 894 }/* End kk */ 895 c_rmtjj[loc] = apindices[jj]; 896 c_rmtaa[loc] = apvaluestmp[jj]; 897 c_rmtc[poj[j]]++; 898 } 899 } 900 flops += voff; 901 } /* End jj */ 902 } /* End j */ 903 } /* End i */ 904 905 ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr); 906 907 ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 908 ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr); 909 910 ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 911 ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 912 ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 913 cd = (Mat_SeqAIJ*)(c->A)->data; 914 co = (Mat_SeqAIJ*)(c->B)->data; 915 916 cmaxr = 0; 917 for (i=0; i<pn; i++) { 918 cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i])); 919 } 920 ierr = PetscCalloc5(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pn,&dcc,pn,&occ);CHKERRQ(ierr); 921 ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 922 ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 923 for (i=0; i<am && pn; i++) { 924 ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 925 nzi = pd->i[i+1] - pd->i[i]; 926 if (!nzi) continue; 927 ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr); 928 voff = 0; 929 ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 930 if (!voff) continue; 931 /* Form C(ii, :) */ 932 pdj = pd->j + pd->i[i]; 933 pda = pd->a + pd->i[i]; 934 for (j=0; j<nzi; j++) { 935 row = pcstart + pdj[j]; 936 for (jj=0; jj<voff; jj++) { 937 apvaluestmp[jj] = apvalues[jj]*pda[j]; 938 } 939 ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr); 940 flops += voff; 941 } 942 } 943 944 ierr = MatGetOwnershipRangeColumn(C,&ccstart,&ccend);CHKERRQ(ierr); 945 ierr = PetscFree5(apindices,apvalues,apvaluestmp,dcc,occ);CHKERRQ(ierr); 946 ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr); 947 ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr); 948 ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 949 ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr); 950 951 /* Add contributions from remote */ 952 for (i = 0; i < pn; i++) { 953 row = i + pcstart; 954 ierr = MatSetValues(C,1,&row,ptap->c_othi[i+1]-ptap->c_othi[i],c_othj+ptap->c_othi[i],c_otha+ptap->c_othi[i],ADD_VALUES);CHKERRQ(ierr); 955 } 956 ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr); 957 958 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 959 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 960 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 961 962 ptap->reuse = MAT_REUSE_MATRIX; 963 964 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 965 if (ptap->freestruct) { 966 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 967 } 968 PetscFunctionReturn(0); 969 } 970 971 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,Mat C) 972 { 973 PetscErrorCode ierr; 974 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 975 Mat_SeqAIJ *cd,*co,*po,*pd; 976 Mat_APMPI *ptap = c->ap; 977 PetscHMapIV hmap; 978 PetscInt i,j,jj,kk,nzi,dnzi,*c_rmtj,voff,*c_othj,pn,pon,pcstart,pcend,row,am,*poj,*pdj,*apindices,cmaxr,*c_rmtc,*c_rmtjj,loc; 979 PetscScalar *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa; 980 PetscLogDouble flops=0.0; 981 MPI_Comm comm; 982 983 PetscFunctionBegin; 984 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 985 if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 986 987 ierr = MatZeroEntries(C);CHKERRQ(ierr); 988 989 /* 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 } 995 996 po = (Mat_SeqAIJ*) p->B->data; 997 pd = (Mat_SeqAIJ*) p->A->data; 998 999 ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 1000 ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 1001 1002 ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr); 1003 ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 1004 ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 1005 cmaxr = 0; 1006 for (i=0; i<pon; i++) { 1007 cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]); 1008 } 1009 cd = (Mat_SeqAIJ*)(c->A)->data; 1010 co = (Mat_SeqAIJ*)(c->B)->data; 1011 for (i=0; i<pn; i++) { 1012 cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i])); 1013 } 1014 ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr); 1015 ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 1016 ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 1017 for (i=0; i<am && (pon || pn); i++) { 1018 ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 1019 nzi = po->i[i+1] - po->i[i]; 1020 dnzi = pd->i[i+1] - pd->i[i]; 1021 if (!nzi && !dnzi) continue; 1022 ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr); 1023 voff = 0; 1024 ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 1025 if (!voff) continue; 1026 1027 /* Form remote C(ii, :) */ 1028 poj = po->j + po->i[i]; 1029 poa = po->a + po->i[i]; 1030 for (j=0; j<nzi; j++) { 1031 c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]]; 1032 c_rmtaa = c_rmta + ptap->c_rmti[poj[j]]; 1033 for (jj=0; jj<voff; jj++) { 1034 apvaluestmp[jj] = apvalues[jj]*poa[j]; 1035 /*If the row is empty */ 1036 if (!c_rmtc[poj[j]]) { 1037 c_rmtjj[jj] = apindices[jj]; 1038 c_rmtaa[jj] = apvaluestmp[jj]; 1039 c_rmtc[poj[j]]++; 1040 } else { 1041 ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr); 1042 if (loc>=0){ /* hit */ 1043 c_rmtaa[loc] += apvaluestmp[jj]; 1044 } else { /* new element */ 1045 loc = -(loc+1); 1046 /* Move data backward */ 1047 for (kk=c_rmtc[poj[j]]; kk>loc; kk--) { 1048 c_rmtjj[kk] = c_rmtjj[kk-1]; 1049 c_rmtaa[kk] = c_rmtaa[kk-1]; 1050 }/* End kk */ 1051 c_rmtjj[loc] = apindices[jj]; 1052 c_rmtaa[loc] = apvaluestmp[jj]; 1053 c_rmtc[poj[j]]++; 1054 } 1055 } 1056 } /* End jj */ 1057 flops += voff; 1058 } /* End j */ 1059 1060 /* Form local C(ii, :) */ 1061 pdj = pd->j + pd->i[i]; 1062 pda = pd->a + pd->i[i]; 1063 for (j=0; j<dnzi; j++) { 1064 row = pcstart + pdj[j]; 1065 for (jj=0; jj<voff; jj++) { 1066 apvaluestmp[jj] = apvalues[jj]*pda[j]; 1067 }/* End kk */ 1068 ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr); 1069 flops += voff; 1070 }/* End j */ 1071 } /* End i */ 1072 1073 ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr); 1074 ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr); 1075 ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr); 1076 1077 ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 1078 ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 1079 ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr); 1080 ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 1081 ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr); 1082 1083 /* Add contributions from remote */ 1084 for (i = 0; i < pn; i++) { 1085 row = i + pcstart; 1086 ierr = MatSetValues(C,1,&row,ptap->c_othi[i+1]-ptap->c_othi[i],c_othj+ptap->c_othi[i],c_otha+ptap->c_othi[i],ADD_VALUES);CHKERRQ(ierr); 1087 } 1088 ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr); 1089 1090 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1091 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1092 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1093 1094 ptap->reuse = MAT_REUSE_MATRIX; 1095 1096 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1097 if (ptap->freestruct) { 1098 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 1099 } 1100 PetscFunctionReturn(0); 1101 } 1102 1103 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,PetscReal fill,Mat *C) 1104 { 1105 Mat_APMPI *ptap; 1106 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1107 MPI_Comm comm; 1108 Mat Cmpi; 1109 Mat_SeqAIJ *pd,*po; 1110 MatType mtype; 1111 PetscSF sf; 1112 PetscSFNode *iremote; 1113 PetscInt rootspacesize,*rootspace,*rootspaceoffsets,nleaves; 1114 const PetscInt *rootdegrees; 1115 PetscHSetI ht,oht,*hta,*hto; 1116 PetscInt pn,pon,*c_rmtc,i,j,nzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets; 1117 PetscInt owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj; 1118 PetscInt nalg=2,alg=0; 1119 PetscBool flg; 1120 const char *algTypes[2] = {"overlapping","merged"}; 1121 PetscErrorCode ierr; 1122 1123 PetscFunctionBegin; 1124 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1125 1126 /* Create symbolic parallel matrix Cmpi */ 1127 ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 1128 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1129 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1130 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1131 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1132 1133 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1134 ptap->reuse = MAT_INITIAL_MATRIX; 1135 ptap->algType = 2; 1136 1137 /* Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1138 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1139 1140 po = (Mat_SeqAIJ*)p->B->data; 1141 pd = (Mat_SeqAIJ*)p->A->data; 1142 1143 /* This equals to the number of offdiag columns in P */ 1144 ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 1145 /* offsets */ 1146 ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr); 1147 /* The number of columns we will send to remote ranks */ 1148 ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr); 1149 ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr); 1150 for (i=0; i<pon; i++) { 1151 ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 1152 } 1153 ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 1154 ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr); 1155 /* Create hash table to merge all columns for C(i, :) */ 1156 ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 1157 1158 ptap->c_rmti[0] = 0; 1159 /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors) */ 1160 for (i=0; i<am && pon; i++) { 1161 /* Form one row of AP */ 1162 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 1163 /* If the off diag is empty, we should not do any calculation */ 1164 nzi = po->i[i+1] - po->i[i]; 1165 if (!nzi) continue; 1166 1167 ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,ht);CHKERRQ(ierr); 1168 ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 1169 /* If AP is empty, just continue */ 1170 if (!htsize) continue; 1171 /* Form C(ii, :) */ 1172 poj = po->j + po->i[i]; 1173 for (j=0; j<nzi; j++) { 1174 ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr); 1175 } 1176 } 1177 1178 for (i=0; i<pon; i++) { 1179 ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 1180 ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize; 1181 c_rmtc[i] = htsize; 1182 } 1183 1184 ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr); 1185 1186 for (i=0; i<pon; i++) { 1187 off = 0; 1188 ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr); 1189 ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 1190 } 1191 ierr = PetscFree(hta);CHKERRQ(ierr); 1192 1193 ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr); 1194 for (i=0; i<pon; i++) { 1195 owner = 0; lidx = 0; 1196 ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 1197 iremote[i].index = lidx; 1198 iremote[i].rank = owner; 1199 } 1200 1201 ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 1202 ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1203 /* Reorder ranks properly so that the data handled by gather and scatter have the same order */ 1204 ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr); 1205 ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 1206 ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1207 /* How many neighbors have contributions to my rows? */ 1208 ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr); 1209 ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr); 1210 rootspacesize = 0; 1211 for (i = 0; i < pn; i++) { 1212 rootspacesize += rootdegrees[i]; 1213 } 1214 ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr); 1215 ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr); 1216 /* Get information from leaves 1217 * Number of columns other people contribute to my rows 1218 * */ 1219 ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 1220 ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 1221 ierr = PetscFree(c_rmtc);CHKERRQ(ierr); 1222 ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr); 1223 /* The number of columns is received for each row */ 1224 ptap->c_othi[0] = 0; 1225 rootspacesize = 0; 1226 rootspaceoffsets[0] = 0; 1227 for (i = 0; i < pn; i++) { 1228 rcvncols = 0; 1229 for (j = 0; j<rootdegrees[i]; j++) { 1230 rcvncols += rootspace[rootspacesize]; 1231 rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize]; 1232 rootspacesize++; 1233 } 1234 ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols; 1235 } 1236 ierr = PetscFree(rootspace);CHKERRQ(ierr); 1237 1238 ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr); 1239 ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 1240 ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 1241 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 1242 ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr); 1243 1244 ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr); 1245 nleaves = 0; 1246 for (i = 0; i<pon; i++) { 1247 owner = 0; lidx = 0; 1248 ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 1249 sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i]; 1250 for (j=0; j<sendncols; j++) { 1251 iremote[nleaves].rank = owner; 1252 iremote[nleaves++].index = c_rmtoffsets[i] + j; 1253 } 1254 } 1255 ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr); 1256 ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr); 1257 1258 ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr); 1259 ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1260 ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr); 1261 /* One to one map */ 1262 ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 1263 1264 ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr); 1265 ierr = PetscHSetICreate(&oht);CHKERRQ(ierr); 1266 ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 1267 ierr = PetscMalloc2(pn,&hta,pn,&hto);CHKERRQ(ierr); 1268 for (i=0; i<pn; i++) { 1269 ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 1270 ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr); 1271 } 1272 /* Work on local part */ 1273 /* 4) Pass 1: Estimate memory for C_loc */ 1274 for (i=0; i<am && pn; i++) { 1275 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 1276 ierr = PetscHSetIClear(oht);CHKERRQ(ierr); 1277 nzi = pd->i[i+1] - pd->i[i]; 1278 if (!nzi) continue; 1279 1280 ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr); 1281 ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 1282 ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr); 1283 if (!(htsize+htosize)) continue; 1284 /* Form C(ii, :) */ 1285 pdj = pd->j + pd->i[i]; 1286 for (j=0; j<nzi; j++) { 1287 ierr = PetscHSetIUpdate(hta[pdj[j]],ht);CHKERRQ(ierr); 1288 ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr); 1289 } 1290 } 1291 1292 ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 1293 ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr); 1294 1295 /* Get remote data */ 1296 ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 1297 ierr = PetscFree(c_rmtj);CHKERRQ(ierr); 1298 1299 for (i = 0; i < pn; i++) { 1300 nzi = ptap->c_othi[i+1] - ptap->c_othi[i]; 1301 rdj = c_othj + ptap->c_othi[i]; 1302 for (j = 0; j < nzi; j++) { 1303 col = rdj[j]; 1304 /* diag part */ 1305 if (col>=pcstart && col<pcend) { 1306 ierr = PetscHSetIAdd(hta[i],col);CHKERRQ(ierr); 1307 } else { /* off diag */ 1308 ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr); 1309 } 1310 } 1311 ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 1312 dnz[i] = htsize; 1313 ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 1314 ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr); 1315 onz[i] = htsize; 1316 ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr); 1317 } 1318 1319 ierr = PetscFree2(hta,hto);CHKERRQ(ierr); 1320 ierr = PetscFree(c_othj);CHKERRQ(ierr); 1321 1322 /* local sizes and preallocation */ 1323 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1324 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 1325 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1326 ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 1327 1328 /* attach the supporting struct to Cmpi for reuse */ 1329 c = (Mat_MPIAIJ*)Cmpi->data; 1330 c->ap = ptap; 1331 ptap->duplicate = Cmpi->ops->duplicate; 1332 ptap->destroy = Cmpi->ops->destroy; 1333 ptap->view = Cmpi->ops->view; 1334 1335 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1336 Cmpi->assembled = PETSC_FALSE; 1337 /* pick an algorithm */ 1338 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 1339 alg = 0; 1340 ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 1341 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1342 switch (alg) { 1343 case 0: 1344 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce; 1345 break; 1346 case 1: 1347 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged; 1348 break; 1349 default: 1350 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n"); 1351 } 1352 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1353 Cmpi->ops->view = MatView_MPIAIJ_PtAP; 1354 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1355 *C = Cmpi; 1356 PetscFunctionReturn(0); 1357 } 1358 1359 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,PetscReal fill,Mat *C) 1360 { 1361 Mat_APMPI *ptap; 1362 Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1363 MPI_Comm comm; 1364 Mat Cmpi; 1365 Mat_SeqAIJ *pd,*po; 1366 MatType mtype; 1367 PetscSF sf; 1368 PetscSFNode *iremote; 1369 PetscInt rootspacesize,*rootspace,*rootspaceoffsets,nleaves; 1370 const PetscInt *rootdegrees; 1371 PetscHSetI ht,oht,*hta,*hto,*htd; 1372 PetscInt pn,pon,*c_rmtc,i,j,nzi,dnzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets; 1373 PetscInt owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj; 1374 PetscInt nalg=2,alg=0; 1375 PetscBool flg; 1376 const char *algTypes[2] = {"merged","overlapping"}; 1377 PetscErrorCode ierr; 1378 1379 PetscFunctionBegin; 1380 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1381 1382 /* Create symbolic parallel matrix Cmpi */ 1383 ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 1384 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1385 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1386 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1387 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1388 1389 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1390 ptap->reuse = MAT_INITIAL_MATRIX; 1391 ptap->algType = 3; 1392 1393 /* 0) Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1394 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1395 1396 po = (Mat_SeqAIJ*)p->B->data; 1397 pd = (Mat_SeqAIJ*)p->A->data; 1398 1399 /* This equals to the number of offdiag columns in P */ 1400 ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 1401 /* offsets */ 1402 ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr); 1403 /* The number of columns we will send to remote ranks */ 1404 ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr); 1405 ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr); 1406 for (i=0; i<pon; i++) { 1407 ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 1408 } 1409 ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 1410 ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr); 1411 /* Create hash table to merge all columns for C(i, :) */ 1412 ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 1413 ierr = PetscHSetICreate(&oht);CHKERRQ(ierr); 1414 ierr = PetscMalloc2(pn,&htd,pn,&hto);CHKERRQ(ierr); 1415 for (i=0; i<pn; i++) { 1416 ierr = PetscHSetICreate(&htd[i]);CHKERRQ(ierr); 1417 ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr); 1418 } 1419 ptap->c_rmti[0] = 0; 1420 /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors) */ 1421 for (i=0; i<am && (pon || pn); i++) { 1422 /* Form one row of AP */ 1423 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 1424 ierr = PetscHSetIClear(oht);CHKERRQ(ierr); 1425 /* If the off diag is empty, we should not do any calculation */ 1426 nzi = po->i[i+1] - po->i[i]; 1427 dnzi = pd->i[i+1] - pd->i[i]; 1428 if (!nzi && !dnzi) continue; 1429 1430 ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr); 1431 ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 1432 ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr); 1433 /* If AP is empty, just continue */ 1434 if (!(htsize+htosize)) continue; 1435 1436 /* Form remote C(ii, :) */ 1437 poj = po->j + po->i[i]; 1438 for (j=0; j<nzi; j++) { 1439 ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr); 1440 ierr = PetscHSetIUpdate(hta[poj[j]],oht);CHKERRQ(ierr); 1441 } 1442 1443 /* Form local C(ii, :) */ 1444 pdj = pd->j + pd->i[i]; 1445 for (j=0; j<dnzi; j++) { 1446 ierr = PetscHSetIUpdate(htd[pdj[j]],ht);CHKERRQ(ierr); 1447 ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr); 1448 } 1449 } 1450 1451 ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 1452 ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr); 1453 1454 for (i=0; i<pon; i++) { 1455 ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 1456 ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize; 1457 c_rmtc[i] = htsize; 1458 } 1459 1460 ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr); 1461 1462 for (i=0; i<pon; i++) { 1463 off = 0; 1464 ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr); 1465 ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 1466 } 1467 ierr = PetscFree(hta);CHKERRQ(ierr); 1468 1469 ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr); 1470 for (i=0; i<pon; i++) { 1471 owner = 0; lidx = 0; 1472 ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 1473 iremote[i].index = lidx; 1474 iremote[i].rank = owner; 1475 } 1476 1477 ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 1478 ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1479 /* Reorder ranks properly so that the data handled by gather and scatter have the same order */ 1480 ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr); 1481 ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 1482 ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 1483 /* How many neighbors have contributions to my rows? */ 1484 ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr); 1485 ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr); 1486 rootspacesize = 0; 1487 for (i = 0; i < pn; i++) { 1488 rootspacesize += rootdegrees[i]; 1489 } 1490 ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr); 1491 ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr); 1492 /* Get information from leaves 1493 * Number of columns other people contribute to my rows 1494 * */ 1495 ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 1496 ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 1497 ierr = PetscFree(c_rmtc);CHKERRQ(ierr); 1498 ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr); 1499 /* The number of columns is received for each row */ 1500 ptap->c_othi[0] = 0; 1501 rootspacesize = 0; 1502 rootspaceoffsets[0] = 0; 1503 for (i = 0; i < pn; i++) { 1504 rcvncols = 0; 1505 for (j = 0; j<rootdegrees[i]; j++) { 1506 rcvncols += rootspace[rootspacesize]; 1507 rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize]; 1508 rootspacesize++; 1509 } 1510 ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols; 1511 } 1512 ierr = PetscFree(rootspace);CHKERRQ(ierr); 1513 1514 ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr); 1515 ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 1516 ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 1517 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 1518 ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr); 1519 1520 ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr); 1521 nleaves = 0; 1522 for (i = 0; i<pon; i++) { 1523 owner = 0; lidx = 0; 1524 ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 1525 sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i]; 1526 for (j=0; j<sendncols; j++) { 1527 iremote[nleaves].rank = owner; 1528 iremote[nleaves++].index = c_rmtoffsets[i] + j; 1529 } 1530 } 1531 ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr); 1532 ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr); 1533 1534 ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr); 1535 ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 1536 ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr); 1537 /* One to one map */ 1538 ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 1539 /* Get remote data */ 1540 ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 1541 ierr = PetscFree(c_rmtj);CHKERRQ(ierr); 1542 ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr); 1543 ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 1544 1545 for (i = 0; i < pn; i++) { 1546 nzi = ptap->c_othi[i+1] - ptap->c_othi[i]; 1547 rdj = c_othj + ptap->c_othi[i]; 1548 for (j = 0; j < nzi; j++) { 1549 col = rdj[j]; 1550 /* diag part */ 1551 if (col>=pcstart && col<pcend) { 1552 ierr = PetscHSetIAdd(htd[i],col);CHKERRQ(ierr); 1553 } else { /* off diag */ 1554 ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr); 1555 } 1556 } 1557 ierr = PetscHSetIGetSize(htd[i],&htsize);CHKERRQ(ierr); 1558 dnz[i] = htsize; 1559 ierr = PetscHSetIDestroy(&htd[i]);CHKERRQ(ierr); 1560 ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr); 1561 onz[i] = htsize; 1562 ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr); 1563 } 1564 1565 ierr = PetscFree2(htd,hto);CHKERRQ(ierr); 1566 ierr = PetscFree(c_othj);CHKERRQ(ierr); 1567 1568 /* local sizes and preallocation */ 1569 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1570 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 1571 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1572 ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 1573 1574 /* attach the supporting struct to Cmpi for reuse */ 1575 c = (Mat_MPIAIJ*)Cmpi->data; 1576 c->ap = ptap; 1577 ptap->duplicate = Cmpi->ops->duplicate; 1578 ptap->destroy = Cmpi->ops->destroy; 1579 ptap->view = Cmpi->ops->view; 1580 1581 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1582 Cmpi->assembled = PETSC_FALSE; 1583 /* pick an algorithm */ 1584 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 1585 alg = 0; 1586 ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 1587 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1588 switch (alg) { 1589 case 0: 1590 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged; 1591 break; 1592 case 1: 1593 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce; 1594 break; 1595 default: 1596 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n"); 1597 } 1598 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1599 Cmpi->ops->view = MatView_MPIAIJ_PtAP; 1600 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1601 *C = Cmpi; 1602 PetscFunctionReturn(0); 1603 } 1604 1605 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 1606 { 1607 PetscErrorCode ierr; 1608 Mat_APMPI *ptap; 1609 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 1610 MPI_Comm comm; 1611 PetscMPIInt size,rank; 1612 Mat Cmpi; 1613 PetscFreeSpaceList free_space=NULL,current_space=NULL; 1614 PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 1615 PetscInt *lnk,i,k,pnz,row,nsend; 1616 PetscBT lnkbt; 1617 PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 1618 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1619 PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 1620 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1621 MPI_Request *swaits,*rwaits; 1622 MPI_Status *sstatus,rstatus; 1623 PetscLayout rowmap; 1624 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 1625 PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 1626 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 1627 Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 1628 PetscScalar *apv; 1629 PetscTable ta; 1630 MatType mtype; 1631 const char *prefix; 1632 #if defined(PETSC_USE_INFO) 1633 PetscReal apfill; 1634 #endif 1635 1636 PetscFunctionBegin; 1637 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1638 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1639 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1640 1641 if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 1642 1643 /* create symbolic parallel matrix Cmpi */ 1644 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1645 ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1646 ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1647 1648 /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 1649 Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 1650 1651 /* create struct Mat_APMPI and attached it to C later */ 1652 ierr = PetscNew(&ptap);CHKERRQ(ierr); 1653 ptap->reuse = MAT_INITIAL_MATRIX; 1654 ptap->algType = 1; 1655 1656 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1657 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1658 /* get P_loc by taking all local rows of P */ 1659 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1660 1661 /* (0) compute Rd = Pd^T, Ro = Po^T */ 1662 /* --------------------------------- */ 1663 ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1664 ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1665 1666 /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 1667 /* ----------------------------------------------------------------- */ 1668 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1669 if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1670 1671 /* create and initialize a linked list */ 1672 ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 1673 MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 1674 MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 1675 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 1676 /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */ 1677 1678 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1679 1680 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 1681 if (ao) { 1682 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 1683 } else { 1684 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 1685 } 1686 current_space = free_space; 1687 nspacedouble = 0; 1688 1689 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 1690 api[0] = 0; 1691 for (i=0; i<am; i++) { 1692 /* diagonal portion: Ad[i,:]*P */ 1693 ai = ad->i; pi = p_loc->i; 1694 nzi = ai[i+1] - ai[i]; 1695 aj = ad->j + ai[i]; 1696 for (j=0; j<nzi; j++) { 1697 row = aj[j]; 1698 pnz = pi[row+1] - pi[row]; 1699 Jptr = p_loc->j + pi[row]; 1700 /* add non-zero cols of P into the sorted linked list lnk */ 1701 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1702 } 1703 /* off-diagonal portion: Ao[i,:]*P */ 1704 if (ao) { 1705 ai = ao->i; pi = p_oth->i; 1706 nzi = ai[i+1] - ai[i]; 1707 aj = ao->j + ai[i]; 1708 for (j=0; j<nzi; j++) { 1709 row = aj[j]; 1710 pnz = pi[row+1] - pi[row]; 1711 Jptr = p_oth->j + pi[row]; 1712 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1713 } 1714 } 1715 apnz = lnk[0]; 1716 api[i+1] = api[i] + apnz; 1717 if (ap_rmax < apnz) ap_rmax = apnz; 1718 1719 /* if free space is not available, double the total space in the list */ 1720 if (current_space->local_remaining<apnz) { 1721 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1722 nspacedouble++; 1723 } 1724 1725 /* Copy data into free space, then initialize lnk */ 1726 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1727 1728 current_space->array += apnz; 1729 current_space->local_used += apnz; 1730 current_space->local_remaining -= apnz; 1731 } 1732 /* Allocate space for apj and apv, initialize apj, and */ 1733 /* destroy list of free space and other temporary array(s) */ 1734 ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 1735 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 1736 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1737 1738 /* Create AP_loc for reuse */ 1739 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 1740 1741 #if defined(PETSC_USE_INFO) 1742 if (ao) { 1743 apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 1744 } else { 1745 apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 1746 } 1747 ptap->AP_loc->info.mallocs = nspacedouble; 1748 ptap->AP_loc->info.fill_ratio_given = fill; 1749 ptap->AP_loc->info.fill_ratio_needed = apfill; 1750 1751 if (api[am]) { 1752 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); 1753 ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 1754 } else { 1755 ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 1756 } 1757 #endif 1758 1759 /* (2-1) compute symbolic Co = Ro*AP_loc */ 1760 /* ------------------------------------ */ 1761 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1762 ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 1763 ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 1764 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 1765 1766 /* (3) send coj of C_oth to other processors */ 1767 /* ------------------------------------------ */ 1768 /* determine row ownership */ 1769 ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 1770 rowmap->n = pn; 1771 rowmap->bs = 1; 1772 ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 1773 owners = rowmap->range; 1774 1775 /* determine the number of messages to send, their lengths */ 1776 ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 1777 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1778 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1779 1780 c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 1781 coi = c_oth->i; coj = c_oth->j; 1782 con = ptap->C_oth->rmap->n; 1783 proc = 0; 1784 for (i=0; i<con; i++) { 1785 while (prmap[i] >= owners[proc+1]) proc++; 1786 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 1787 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 1788 } 1789 1790 len = 0; /* max length of buf_si[], see (4) */ 1791 owners_co[0] = 0; 1792 nsend = 0; 1793 for (proc=0; proc<size; proc++) { 1794 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1795 if (len_s[proc]) { 1796 nsend++; 1797 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 1798 len += len_si[proc]; 1799 } 1800 } 1801 1802 /* determine the number and length of messages to receive for coi and coj */ 1803 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 1804 ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 1805 1806 /* post the Irecv and Isend of coj */ 1807 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1808 ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1809 ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 1810 for (proc=0, k=0; proc<size; proc++) { 1811 if (!len_s[proc]) continue; 1812 i = owners_co[proc]; 1813 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1814 k++; 1815 } 1816 1817 /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 1818 /* ---------------------------------------- */ 1819 ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 1820 ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 1821 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 1822 c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 1823 1824 /* receives coj are complete */ 1825 for (i=0; i<nrecv; i++) { 1826 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1827 } 1828 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1829 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1830 1831 /* add received column indices into ta to update Crmax */ 1832 for (k=0; k<nrecv; k++) {/* k-th received message */ 1833 Jptr = buf_rj[k]; 1834 for (j=0; j<len_r[k]; j++) { 1835 ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1836 } 1837 } 1838 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 1839 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 1840 1841 /* (4) send and recv coi */ 1842 /*-----------------------*/ 1843 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1844 ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1845 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1846 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1847 for (proc=0,k=0; proc<size; proc++) { 1848 if (!len_s[proc]) continue; 1849 /* form outgoing message for i-structure: 1850 buf_si[0]: nrows to be sent 1851 [1:nrows]: row index (global) 1852 [nrows+1:2*nrows+1]: i-structure index 1853 */ 1854 /*-------------------------------------------*/ 1855 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 1856 buf_si_i = buf_si + nrows+1; 1857 buf_si[0] = nrows; 1858 buf_si_i[0] = 0; 1859 nrows = 0; 1860 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1861 nzi = coi[i+1] - coi[i]; 1862 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1863 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1864 nrows++; 1865 } 1866 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1867 k++; 1868 buf_si += len_si[proc]; 1869 } 1870 for (i=0; i<nrecv; i++) { 1871 ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1872 } 1873 ierr = PetscFree(rwaits);CHKERRQ(ierr); 1874 if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1875 1876 ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 1877 ierr = PetscFree(len_ri);CHKERRQ(ierr); 1878 ierr = PetscFree(swaits);CHKERRQ(ierr); 1879 ierr = PetscFree(buf_s);CHKERRQ(ierr); 1880 1881 /* (5) compute the local portion of Cmpi */ 1882 /* ------------------------------------------ */ 1883 /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 1884 ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 1885 current_space = free_space; 1886 1887 ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 1888 for (k=0; k<nrecv; k++) { 1889 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1890 nrows = *buf_ri_k[k]; 1891 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1892 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1893 } 1894 1895 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 1896 ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1897 for (i=0; i<pn; i++) { 1898 /* add C_loc into Cmpi */ 1899 nzi = c_loc->i[i+1] - c_loc->i[i]; 1900 Jptr = c_loc->j + c_loc->i[i]; 1901 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1902 1903 /* add received col data into lnk */ 1904 for (k=0; k<nrecv; k++) { /* k-th received message */ 1905 if (i == *nextrow[k]) { /* i-th row */ 1906 nzi = *(nextci[k]+1) - *nextci[k]; 1907 Jptr = buf_rj[k] + *nextci[k]; 1908 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1909 nextrow[k]++; nextci[k]++; 1910 } 1911 } 1912 nzi = lnk[0]; 1913 1914 /* copy data into free space, then initialize lnk */ 1915 ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1916 ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 1917 } 1918 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1919 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1920 ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 1921 1922 /* local sizes and preallocation */ 1923 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1924 if (P->cmap->bs > 0) { 1925 ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr); 1926 ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr); 1927 } 1928 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1929 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1930 1931 /* members in merge */ 1932 ierr = PetscFree(id_r);CHKERRQ(ierr); 1933 ierr = PetscFree(len_r);CHKERRQ(ierr); 1934 ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 1935 ierr = PetscFree(buf_ri);CHKERRQ(ierr); 1936 ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 1937 ierr = PetscFree(buf_rj);CHKERRQ(ierr); 1938 ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 1939 1940 /* attach the supporting struct to Cmpi for reuse */ 1941 c = (Mat_MPIAIJ*)Cmpi->data; 1942 c->ap = ptap; 1943 ptap->duplicate = Cmpi->ops->duplicate; 1944 ptap->destroy = Cmpi->ops->destroy; 1945 ptap->view = Cmpi->ops->view; 1946 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 1947 1948 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1949 Cmpi->assembled = PETSC_FALSE; 1950 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1951 Cmpi->ops->view = MatView_MPIAIJ_PtAP; 1952 Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1953 *C = Cmpi; 1954 PetscFunctionReturn(0); 1955 } 1956 1957 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 1958 { 1959 PetscErrorCode ierr; 1960 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1961 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1962 Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 1963 Mat_APMPI *ptap = c->ap; 1964 Mat AP_loc,C_loc,C_oth; 1965 PetscInt i,rstart,rend,cm,ncols,row; 1966 PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 1967 PetscScalar *apa; 1968 const PetscInt *cols; 1969 const PetscScalar *vals; 1970 1971 PetscFunctionBegin; 1972 if (!ptap->AP_loc) { 1973 MPI_Comm comm; 1974 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1975 SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1976 } 1977 1978 ierr = MatZeroEntries(C);CHKERRQ(ierr); 1979 /* 1) get R = Pd^T,Ro = Po^T */ 1980 if (ptap->reuse == MAT_REUSE_MATRIX) { 1981 ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1982 ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1983 } 1984 1985 /* 2) get AP_loc */ 1986 AP_loc = ptap->AP_loc; 1987 ap = (Mat_SeqAIJ*)AP_loc->data; 1988 1989 /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 1990 /*-----------------------------------------------------*/ 1991 if (ptap->reuse == MAT_REUSE_MATRIX) { 1992 /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 1993 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1994 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1995 } 1996 1997 /* 2-2) compute numeric A_loc*P - dominating part */ 1998 /* ---------------------------------------------- */ 1999 /* get data from symbolic products */ 2000 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 2001 if (ptap->P_oth) { 2002 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 2003 } 2004 apa = ptap->apa; 2005 api = ap->i; 2006 apj = ap->j; 2007 for (i=0; i<am; i++) { 2008 /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 2009 AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 2010 apnz = api[i+1] - api[i]; 2011 for (j=0; j<apnz; j++) { 2012 col = apj[j+api[i]]; 2013 ap->a[j+ap->i[i]] = apa[col]; 2014 apa[col] = 0.0; 2015 } 2016 } 2017 2018 /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 2019 ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 2020 ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 2021 C_loc = ptap->C_loc; 2022 C_oth = ptap->C_oth; 2023 2024 /* add C_loc and Co to to C */ 2025 ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 2026 2027 /* C_loc -> C */ 2028 cm = C_loc->rmap->N; 2029 c_seq = (Mat_SeqAIJ*)C_loc->data; 2030 cols = c_seq->j; 2031 vals = c_seq->a; 2032 2033 2034 /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 2035 /* when there are no off-processor parts. */ 2036 /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 2037 /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 2038 /* a table, and other, more complex stuff has to be done. */ 2039 if (C->assembled) { 2040 C->was_assembled = PETSC_TRUE; 2041 C->assembled = PETSC_FALSE; 2042 } 2043 if (C->was_assembled) { 2044 for (i=0; i<cm; i++) { 2045 ncols = c_seq->i[i+1] - c_seq->i[i]; 2046 row = rstart + i; 2047 ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 2048 cols += ncols; vals += ncols; 2049 } 2050 } else { 2051 ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 2052 } 2053 2054 /* Co -> C, off-processor part */ 2055 cm = C_oth->rmap->N; 2056 c_seq = (Mat_SeqAIJ*)C_oth->data; 2057 cols = c_seq->j; 2058 vals = c_seq->a; 2059 for (i=0; i<cm; i++) { 2060 ncols = c_seq->i[i+1] - c_seq->i[i]; 2061 row = p->garray[i]; 2062 ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 2063 cols += ncols; vals += ncols; 2064 } 2065 2066 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2067 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2068 2069 ptap->reuse = MAT_REUSE_MATRIX; 2070 2071 /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 2072 if (ptap->freestruct) { 2073 ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 2074 } 2075 PetscFunctionReturn(0); 2076 } 2077