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