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 12 EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP" 15 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 16 { 17 PetscErrorCode ierr; 18 Mat_MatMatMultMPI *mult; 19 PetscObjectContainer container; 20 21 PetscFunctionBegin; 22 ierr = PetscObjectQuery((PetscObject)A,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); 23 if (container) { 24 ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 25 ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr); 26 ierr = MatDestroy(mult->B_oth);CHKERRQ(ierr); 27 if (mult->abi) ierr = PetscFree(mult->abi);CHKERRQ(ierr); 28 if (mult->abj) ierr = PetscFree(mult->abj);CHKERRQ(ierr); 29 if (mult->startsj) ierr = PetscFree(mult->startsj);CHKERRQ(ierr); 30 if (mult->bufa) ierr = PetscFree(mult->bufa);CHKERRQ(ierr); 31 ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 32 ierr = PetscObjectCompose((PetscObject)A,"MatPtAPMPI",0);CHKERRQ(ierr); 33 } 34 ierr = PetscFree(mult);CHKERRQ(ierr); 35 36 ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 37 PetscFunctionReturn(0); 38 } 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 42 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 43 { 44 PetscErrorCode ierr; 45 Mat_MatMatMultMPI *mult; 46 PetscObjectContainer container; 47 48 PetscFunctionBegin; 49 if (scall == MAT_INITIAL_MATRIX){ 50 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 51 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 52 mult->B_loc=PETSC_NULL; mult->B_oth=PETSC_NULL; 53 mult->abi=PETSC_NULL; mult->abj=PETSC_NULL; 54 mult->abnz_max = 0; /* symbolic A*P is not done yet */ 55 56 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 57 ierr = MatGetBrowsOfAoCols(A,P,scall,&mult->startsj,&mult->bufa,&mult->B_oth);CHKERRQ(ierr); 58 59 /* get P_loc by taking all local rows of P */ 60 ierr = MatGetLocalMat(P,scall,&mult->B_loc);CHKERRQ(ierr); 61 62 /* attach the supporting struct to P for reuse */ 63 P->ops->destroy = MatDestroy_MPIAIJ_PtAP; 64 ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 65 ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr); 66 ierr = PetscObjectCompose((PetscObject)P,"MatPtAPMPI",(PetscObject)container);CHKERRQ(ierr); 67 68 /* now, compute symbolic local P^T*A*P */ 69 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */ 70 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 71 } else if (scall == MAT_REUSE_MATRIX){ 72 ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); 73 if (container) { 74 ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 75 } else { 76 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 77 } 78 79 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 80 ierr = MatGetBrowsOfAoCols(A,P,scall,&mult->startsj,&mult->bufa,&mult->B_oth);CHKERRQ(ierr); 81 82 /* get P_loc by taking all local rows of P */ 83 ierr = MatGetLocalMat(P,scall,&mult->B_loc);CHKERRQ(ierr); 84 85 } else { 86 SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall); 87 } 88 /* now, compute numeric local P^T*A*P */ 89 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 90 ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 91 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 92 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 98 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 99 { 100 PetscErrorCode ierr; 101 Mat P_loc,P_oth,B_mpi; 102 Mat_MatMatMultMPI *ap; 103 PetscObjectContainer container; 104 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 105 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 106 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 107 Mat_SeqAIJ *p_loc,*p_oth; 108 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 109 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz; 110 PetscInt nlnk,*lnk,*owners_co,*coi,*coj,i,j,k,pnz,row; 111 PetscInt am=A->m,pN=P->N,pn=P->n; 112 PetscBT lnkbt; 113 MPI_Comm comm=A->comm; 114 PetscMPIInt size,rank,tag,*len_si,*len_s,*len_ri; 115 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 116 PetscInt len,proc,*dnz,*onz,*owners; 117 PetscInt nzi,*bi,*bj; 118 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 119 MPI_Request *swaits,*rwaits; 120 MPI_Status *sstatus,rstatus; 121 Mat_Merge_SeqsToMPI *merge; 122 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon; 123 124 PetscFunctionBegin; 125 ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); 126 if (container) { 127 ierr = PetscObjectContainerGetPointer(container,(void **)&ap);CHKERRQ(ierr); 128 } else { 129 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 130 } 131 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 132 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 133 134 /* get data from P_loc and P_oth */ 135 P_loc=ap->B_loc; P_oth=ap->B_oth; 136 p_loc=(Mat_SeqAIJ*)P_loc->data; p_oth=(Mat_SeqAIJ*)P_oth->data; 137 pi_loc=p_loc->i; pj_loc=p_loc->j; pi_oth=p_oth->i; pj_oth=p_oth->j; 138 139 /* first, compute symbolic AP = A_loc*P */ 140 /*--------------------------------------*/ 141 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 142 ap->abi = api; 143 api[0] = 0; 144 145 /* create and initialize a linked list */ 146 nlnk = pN+1; 147 ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 148 149 /* Initial FreeSpace size is fill*nnz(A) */ 150 ierr = GetMoreSpace((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr); 151 current_space = free_space; 152 153 for (i=0;i<am;i++) { 154 apnz = 0; 155 /* diagonal portion of A */ 156 nzi = adi[i+1] - adi[i]; 157 for (j=0; j<nzi; j++){ 158 row = *adj++; 159 pnz = pi_loc[row+1] - pi_loc[row]; 160 Jptr = pj_loc + pi_loc[row]; 161 /* add non-zero cols of P into the sorted linked list lnk */ 162 ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 163 apnz += nlnk; 164 } 165 /* off-diagonal portion of A */ 166 nzi = aoi[i+1] - aoi[i]; 167 for (j=0; j<nzi; j++){ 168 row = *aoj++; 169 pnz = pi_oth[row+1] - pi_oth[row]; 170 Jptr = pj_oth + pi_oth[row]; 171 ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 172 apnz += nlnk; 173 } 174 175 api[i+1] = api[i] + apnz; 176 if (ap->abnz_max < apnz) ap->abnz_max = apnz; 177 178 /* if free space is not available, double the total space in the list */ 179 if (current_space->local_remaining<apnz) { 180 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 181 } 182 183 /* Copy data into free space, then initialize lnk */ 184 ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 185 current_space->array += apnz; 186 current_space->local_used += apnz; 187 current_space->local_remaining -= apnz; 188 } 189 /* Allocate space for apj, initialize apj, and */ 190 /* destroy list of free space and other temporary array(s) */ 191 ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr); 192 apj = ap->abj; 193 ierr = MakeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr); 194 195 /* determine symbolic Co=(p->B)^T*AP - send to others */ 196 /*----------------------------------------------------*/ 197 ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 198 199 /* then, compute symbolic Co = (p->B)^T*AP */ 200 pon = (p->B)->n; /* total num of rows to be sent to other processors 201 >= (num of nonzero rows of C_seq) - pn */ 202 ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 203 coi[0] = 0; 204 205 /* set initial free space to be 3*pon*max( nnz(AP) per row) */ 206 nnz = 3*pon*ap->abnz_max + 1; 207 ierr = GetMoreSpace(nnz,&free_space); 208 current_space = free_space; 209 210 for (i=0; i<pon; i++) { 211 nnz = 0; 212 pnz = poti[i+1] - poti[i]; 213 j = pnz; 214 ptJ = potj + poti[i+1]; 215 while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 216 j--; ptJ--; 217 row = *ptJ; /* row of AP == col of Pot */ 218 apnz = api[row+1] - api[row]; 219 Jptr = apj + api[row]; 220 /* add non-zero cols of AP into the sorted linked list lnk */ 221 ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 222 nnz += nlnk; 223 } 224 225 /* If free space is not available, double the total space in the list */ 226 if (current_space->local_remaining<nnz) { 227 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 228 } 229 230 /* Copy data into free space, and zero out denserows */ 231 ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 232 current_space->array += nnz; 233 current_space->local_used += nnz; 234 current_space->local_remaining -= nnz; 235 coi[i+1] = coi[i] + nnz; 236 } 237 ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 238 ierr = MakeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 239 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 240 241 /* send j-array (coj) of Co to other processors */ 242 /*----------------------------------------------*/ 243 /* determine row ownership */ 244 ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 245 ierr = PetscMapCreate(comm,&merge->rowmap);CHKERRQ(ierr); 246 ierr = PetscMapSetLocalSize(merge->rowmap,pn);CHKERRQ(ierr); 247 ierr = PetscMapSetType(merge->rowmap,MAP_MPI);CHKERRQ(ierr); 248 ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 249 250 /* determine the number of messages to send, their lengths */ 251 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 252 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 253 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 254 len_s = merge->len_s; 255 merge->nsend = 0; 256 257 ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 258 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 259 260 proc = 0; 261 for (i=0; i<pon; i++){ 262 while (prmap[i] >= owners[proc+1]) proc++; 263 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 264 len_s[proc] += coi[i+1] - coi[i]; 265 } 266 267 len = 0; /* max length of buf_si[] */ 268 owners_co[0] = 0; 269 for (proc=0; proc<size; proc++){ 270 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 271 if (len_si[proc]){ 272 merge->nsend++; 273 len_si[proc] = 2*(len_si[proc] + 1); 274 len += len_si[proc]; 275 } 276 } 277 278 /* determine the number and length of messages to receive for coi and coj */ 279 ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 280 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 281 282 /* post the Irecv and Isend of coj */ 283 ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tag);CHKERRQ(ierr); 284 ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 285 286 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 287 288 for (proc=0, k=0; proc<size; proc++){ 289 if (!len_s[proc]) continue; 290 i = owners_co[proc]; 291 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 292 k++; 293 } 294 295 /* receives and sends of coj are complete */ 296 ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 297 i = merge->nrecv; 298 while (i--) { 299 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 300 } 301 ierr = PetscFree(rwaits);CHKERRQ(ierr); 302 if (merge->nsend){ 303 ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 304 } 305 306 /* send and recv coi */ 307 /*-------------------*/ 308 ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 309 310 ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 311 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 312 for (proc=0,k=0; proc<size; proc++){ 313 if (!len_s[proc]) continue; 314 /* form outgoing message for i-structure: 315 buf_si[0]: nrows to be sent 316 [1:nrows]: row index (global) 317 [nrows+1:2*nrows+1]: i-structure index 318 */ 319 /*-------------------------------------------*/ 320 nrows = len_si[proc]/2 - 1; 321 buf_si_i = buf_si + nrows+1; 322 buf_si[0] = nrows; 323 buf_si_i[0] = 0; 324 nrows = 0; 325 for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 326 nzi = coi[i+1] - coi[i]; 327 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 328 buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 329 nrows++; 330 } 331 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 332 k++; 333 buf_si += len_si[proc]; 334 } 335 i = merge->nrecv; 336 while (i--) { 337 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 338 } 339 ierr = PetscFree(rwaits);CHKERRQ(ierr); 340 if (merge->nsend){ 341 ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 342 } 343 344 ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr); 345 for (i=0; i<merge->nrecv; i++){ 346 ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);CHKERRQ(ierr); 347 } 348 349 ierr = PetscFree(len_si);CHKERRQ(ierr); 350 ierr = PetscFree(len_ri);CHKERRQ(ierr); 351 ierr = PetscFree(swaits);CHKERRQ(ierr); 352 ierr = PetscFree(sstatus);CHKERRQ(ierr); 353 ierr = PetscFree(buf_s);CHKERRQ(ierr); 354 355 /* compute the local portion of C (mpi mat) */ 356 /*------------------------------------------*/ 357 ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 358 359 /* allocate bi array and free space for accumulating nonzero column info */ 360 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 361 bi[0] = 0; 362 363 /* set initial free space to be 3*pn*max( nnz(AP) per row) */ 364 nnz = 3*pn*ap->abnz_max + 1; 365 ierr = GetMoreSpace(nnz,&free_space); 366 current_space = free_space; 367 368 ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 369 nextrow = buf_ri_k + merge->nrecv; 370 nextci = nextrow + merge->nrecv; 371 for (k=0; k<merge->nrecv; k++){ 372 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 373 nrows = *buf_ri_k[k]; 374 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 375 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 376 } 377 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 378 for (i=0; i<pn; i++) { 379 /* add pdt[i,:]*AP into lnk */ 380 nnz = 0; 381 pnz = pdti[i+1] - pdti[i]; 382 j = pnz; 383 ptJ = pdtj + pdti[i+1]; 384 while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 385 j--; ptJ--; 386 row = *ptJ; /* row of AP == col of Pt */ 387 apnz = api[row+1] - api[row]; 388 Jptr = apj + api[row]; 389 /* add non-zero cols of AP into the sorted linked list lnk */ 390 ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 391 nnz += nlnk; 392 } 393 /* add received col data into lnk */ 394 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 395 if (i == *nextrow[k]) { /* i-th row */ 396 nzi = *(nextci[k]+1) - *nextci[k]; 397 Jptr = buf_rj[k] + *nextci[k]; 398 ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 399 nnz += nlnk; 400 nextrow[k]++; nextci[k]++; 401 } 402 } 403 404 /* if free space is not available, make more free space */ 405 if (current_space->local_remaining<nnz) { 406 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 407 } 408 /* copy data into free space, then initialize lnk */ 409 ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 410 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 411 current_space->array += nnz; 412 current_space->local_used += nnz; 413 current_space->local_remaining -= nnz; 414 bi[i+1] = bi[i] + nnz; 415 } 416 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 417 ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 418 419 ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 420 ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 421 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 422 423 /* create symbolic parallel matrix B_mpi */ 424 /*---------------------------------------*/ 425 ierr = MatCreate(comm,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE,&B_mpi);CHKERRQ(ierr); 426 ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr); 427 ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr); 428 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 429 430 /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 431 B_mpi->assembled = PETSC_FALSE; 432 B_mpi->ops->destroy = MatDestroy_MPIAIJ_SeqsToMPI; 433 merge->bi = bi; 434 merge->bj = bj; 435 merge->coi = coi; 436 merge->coj = coj; 437 merge->buf_ri = buf_ri; 438 merge->buf_rj = buf_rj; 439 merge->owners_co = owners_co; 440 441 /* attach the supporting struct to B_mpi for reuse */ 442 ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 443 ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr); 444 ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr); 445 *C = B_mpi; 446 447 PetscFunctionReturn(0); 448 } 449 450 #undef __FUNCT__ 451 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 452 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 453 { 454 PetscErrorCode ierr; 455 Mat_Merge_SeqsToMPI *merge; 456 Mat_MatMatMultMPI *ap; 457 PetscObjectContainer cont_merge,cont_ptap; 458 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 459 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 460 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 461 Mat_SeqAIJ *p_loc,*p_oth; 462 PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0; 463 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 464 PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 465 MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth; 466 PetscInt am=A->m,cm=C->m,pon=(p->B)->n; 467 MPI_Comm comm=C->comm; 468 PetscMPIInt size,rank,taga,*len_s; 469 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 470 PetscInt **buf_ri,**buf_rj; 471 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 472 MPI_Request *s_waits,*r_waits; 473 MPI_Status *status; 474 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 475 PetscInt *api,*apj,*coi,*coj; 476 PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=p->cstart,pcend=p->cend; 477 478 PetscFunctionBegin; 479 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 480 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 481 482 ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr); 483 if (cont_merge) { 484 ierr = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr); 485 } else { 486 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container"); 487 } 488 ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr); 489 if (cont_ptap) { 490 ierr = PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr); 491 } else { 492 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 493 } 494 495 /* get data from symbolic products */ 496 p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 497 p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 498 pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc; 499 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 500 501 coi = merge->coi; coj = merge->coj; 502 ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 503 ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 504 505 bi = merge->bi; bj = merge->bj; 506 ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 507 ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 508 ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 509 510 /* get data from symbolic A*P */ 511 ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr); 512 ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr); 513 514 /* compute numeric C_seq=P_loc^T*A_loc*P */ 515 api = ap->abi; apj = ap->abj; 516 for (i=0;i<am;i++) { 517 /* form i-th sparse row of A*P */ 518 apnz = api[i+1] - api[i]; 519 apJ = apj + api[i]; 520 /* diagonal portion of A */ 521 anz = adi[i+1] - adi[i]; 522 for (j=0;j<anz;j++) { 523 row = *adj++; 524 pnz = pi_loc[row+1] - pi_loc[row]; 525 pj = pj_loc + pi_loc[row]; 526 pa = pa_loc + pi_loc[row]; 527 nextp = 0; 528 for (k=0; nextp<pnz; k++) { 529 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 530 apa[k] += (*ada)*pa[nextp++]; 531 } 532 } 533 flops += 2*pnz; 534 ada++; 535 } 536 /* off-diagonal portion of A */ 537 anz = aoi[i+1] - aoi[i]; 538 for (j=0; j<anz; j++) { 539 row = *aoj++; 540 pnz = pi_oth[row+1] - pi_oth[row]; 541 pj = pj_oth + pi_oth[row]; 542 pa = pa_oth + pi_oth[row]; 543 nextp = 0; 544 for (k=0; nextp<pnz; k++) { 545 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 546 apa[k] += (*aoa)*pa[nextp++]; 547 } 548 } 549 flops += 2*pnz; 550 aoa++; 551 } 552 553 /* Compute P_loc[i,:]^T*AP[i,:] using outer product */ 554 pnz = pi_loc[i+1] - pi_loc[i]; 555 for (j=0; j<pnz; j++) { 556 nextap = 0; 557 row = *pJ++; /* global index */ 558 if (row < pcstart || row >=pcend) { /* put the value into Co */ 559 cj = coj + coi[*poJ]; 560 ca = coa + coi[*poJ++]; 561 } else { /* put the value into Cd */ 562 cj = bj + bi[*pdJ]; 563 ca = ba + bi[*pdJ++]; 564 } 565 for (k=0; nextap<apnz; k++) { 566 if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++]; 567 } 568 flops += 2*apnz; 569 pA++; 570 } 571 572 /* zero the current row info for A*P */ 573 ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 574 } 575 ierr = PetscFree(apa);CHKERRQ(ierr); 576 577 /* send and recv matrix values */ 578 /*-----------------------------*/ 579 buf_ri = merge->buf_ri; 580 buf_rj = merge->buf_rj; 581 len_s = merge->len_s; 582 ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&taga);CHKERRQ(ierr); 583 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 584 585 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr); 586 for (proc=0,k=0; proc<size; proc++){ 587 if (!len_s[proc]) continue; 588 i = merge->owners_co[proc]; 589 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 590 k++; 591 } 592 ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr); 593 ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr); 594 ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr); 595 ierr = PetscFree(status);CHKERRQ(ierr); 596 597 ierr = PetscFree(s_waits);CHKERRQ(ierr); 598 ierr = PetscFree(r_waits);CHKERRQ(ierr); 599 ierr = PetscFree(coa);CHKERRQ(ierr); 600 601 /* insert local and received values into C */ 602 /*-----------------------------------------*/ 603 ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 604 nextrow = buf_ri_k + merge->nrecv; 605 nextci = nextrow + merge->nrecv; 606 607 for (k=0; k<merge->nrecv; k++){ 608 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 609 nrows = *(buf_ri_k[k]); 610 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 611 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 612 } 613 614 for (i=0; i<cm; i++) { 615 row = owners[rank] + i; /* global row index of C_seq */ 616 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 617 ba_i = ba + bi[i]; 618 bnz = bi[i+1] - bi[i]; 619 /* add received vals into ba */ 620 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 621 /* i-th row */ 622 if (i == *nextrow[k]) { 623 cnz = *(nextci[k]+1) - *nextci[k]; 624 cj = buf_rj[k] + *(nextci[k]); 625 ca = abuf_r[k] + *(nextci[k]); 626 nextcj = 0; 627 for (j=0; nextcj<cnz; j++){ 628 if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 629 ba_i[j] += ca[nextcj++]; 630 } 631 } 632 nextrow[k]++; nextci[k]++; 633 } 634 } 635 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 636 flops += 2*cnz; 637 } 638 ierr = MatSetOption(C,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* -cause delay? */ 639 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 640 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 641 642 ierr = PetscFree(ba);CHKERRQ(ierr); 643 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 644 ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 645 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 646 647 PetscFunctionReturn(0); 648 } 649