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