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