1 /*$Id: mmdense.c,v 1.28 2000/04/09 04:35:58 bsmith Exp bsmith $*/ 2 3 /* 4 Support for the parallel dense matrix vector multiply 5 */ 6 #include "src/mat/impls/dense/mpi/mpidense.h" 7 #include "src/vec/vecimpl.h" 8 9 #undef __FUNC__ 10 #define __FUNC__ /*<a name=""></a>*/"MatSetUpMultiply_MPIDense" 11 int MatSetUpMultiply_MPIDense(Mat mat) 12 { 13 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 14 int ierr; 15 IS tofrom; 16 Vec gvec; 17 18 PetscFunctionBegin; 19 /* Create local vector that is used to scatter into */ 20 ierr = VecCreateSeq(PETSC_COMM_SELF,mdn->N,&mdn->lvec);CHKERRQ(ierr); 21 22 /* Create temporary index set for building scatter gather */ 23 ierr = ISCreateStride(PETSC_COMM_SELF,mdn->N,0,1,&tofrom);CHKERRQ(ierr); 24 25 /* Create temporary global vector to generate scatter context */ 26 /* n = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank]; */ 27 28 ierr = VecCreateMPI(mat->comm,mdn->nvec,mdn->N,&gvec);CHKERRQ(ierr); 29 30 /* Generate the scatter context */ 31 ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx);CHKERRQ(ierr); 32 PLogObjectParent(mat,mdn->Mvctx); 33 PLogObjectParent(mat,mdn->lvec); 34 PLogObjectParent(mat,tofrom); 35 PLogObjectParent(mat,gvec); 36 37 ierr = ISDestroy(tofrom);CHKERRQ(ierr); 38 ierr = VecDestroy(gvec);CHKERRQ(ierr); 39 PetscFunctionReturn(0); 40 } 41 42 extern int MatGetSubMatrices_MPIDense_Local(Mat,int,IS*,IS*,MatReuse,Mat*); 43 #undef __FUNC__ 44 #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrices_MPIDense" 45 int MatGetSubMatrices_MPIDense(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat **submat) 46 { 47 Mat_MPIDense *c = (Mat_MPIDense*)C->data; 48 int nmax,nstages_local,nstages,i,pos,max_no,ierr; 49 50 PetscFunctionBegin; 51 /* Allocate memory to hold all the submatrices */ 52 if (scall != MAT_REUSE_MATRIX) { 53 *submat = (Mat *)PetscMalloc((ismax+1)*sizeof(Mat));CHKPTRQ(*submat); 54 } 55 /* Determine the number of stages through which submatrices are done */ 56 nmax = 20*1000000 / (c->N * sizeof(int)); 57 if (!nmax) nmax = 1; 58 nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 59 60 /* Make sure every processor loops through the nstages */ 61 ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr); 62 63 64 for (i=0,pos=0; i<nstages; i++) { 65 if (pos+nmax <= ismax) max_no = nmax; 66 else if (pos == ismax) max_no = 0; 67 else max_no = ismax-pos; 68 ierr = MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr); 69 pos += max_no; 70 } 71 PetscFunctionReturn(0); 72 } 73 /* -------------------------------------------------------------------------*/ 74 #undef __FUNC__ 75 #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrices_MPIDense_Local" 76 int MatGetSubMatrices_MPIDense_Local(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat *submats) 77 { 78 Mat_MPIDense *c = (Mat_MPIDense*)C->data; 79 Mat A = c->A; 80 Mat_SeqDense *a = (Mat_SeqDense*)A->data,*mat; 81 int N = c->N,rstart = c->rstart,count; 82 int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size; 83 int **sbuf1,rank,m,i,j,k,l,ct1,ierr,**rbuf1,row,proc; 84 int nrqs,msz,**ptr,index,*ctr,*pa,*tmp,bsz,nrqr; 85 int is_no,jmax,*irow_i,**rmap,*rmap_i; 86 int len,ctr_j,*sbuf1_j,*rbuf1_i; 87 int tag0,tag1; 88 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 89 MPI_Status *r_status1,*r_status2,*s_status1,*s_status2; 90 MPI_Comm comm; 91 Scalar **rbuf2,**sbuf2; 92 93 PetscFunctionBegin; 94 comm = C->comm; 95 tag0 = C->tag; 96 size = c->size; 97 rank = c->rank; 98 m = c->M; 99 100 /* Get some new tags to keep the communication clean */ 101 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 102 103 /* Check if the col indices are sorted */ 104 for (i=0; i<ismax; i++) { 105 ierr = ISSorted(isrow[i],(PetscTruth*)&j);CHKERRQ(ierr); 106 if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted"); 107 ierr = ISSorted(iscol[i],(PetscTruth*)&j);CHKERRQ(ierr); 108 if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted"); 109 } 110 111 len = 2*ismax*(sizeof(int *)+sizeof(int)) + (m+1)*sizeof(int); 112 irow = (int **)PetscMalloc(len);CHKPTRQ(irow); 113 icol = irow + ismax; 114 nrow = (int*)(icol + ismax); 115 ncol = nrow + ismax; 116 rtable = ncol + ismax; 117 118 for (i=0; i<ismax; i++) { 119 ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 120 ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 121 ierr = ISGetSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 122 ierr = ISGetSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 123 } 124 125 /* Create hash table for the mapping :row -> proc*/ 126 for (i=0,j=0; i<size; i++) { 127 jmax = c->rowners[i+1]; 128 for (; j<jmax; j++) { 129 rtable[j] = i; 130 } 131 } 132 133 /* evaluate communication - mesg to who,length of mesg, and buffer space 134 required. Based on this, buffers are allocated, and data copied into them*/ 135 w1 = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */ 136 w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 137 w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 138 w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 139 ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/ 140 for (i=0; i<ismax; i++) { 141 ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/ 142 jmax = nrow[i]; 143 irow_i = irow[i]; 144 for (j=0; j<jmax; j++) { 145 row = irow_i[j]; 146 proc = rtable[row]; 147 w4[proc]++; 148 } 149 for (j=0; j<size; j++) { 150 if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 151 } 152 } 153 154 nrqs = 0; /* no of outgoing messages */ 155 msz = 0; /* total mesg length (for all procs) */ 156 w1[rank] = 0; /* no mesg sent to self */ 157 w3[rank] = 0; 158 for (i=0; i<size; i++) { 159 if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 160 } 161 pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/ 162 for (i=0,j=0; i<size; i++) { 163 if (w1[i]) { pa[j] = i; j++; } 164 } 165 166 /* Each message would have a header = 1 + 2*(no of IS) + data */ 167 for (i=0; i<nrqs; i++) { 168 j = pa[i]; 169 w1[j] += w2[j] + 2* w3[j]; 170 msz += w1[j]; 171 } 172 /* Do a global reduction to determine how many messages to expect*/ 173 { 174 int *rw1; 175 rw1 = (int *)PetscMalloc(2*size*sizeof(int));CHKPTRQ(rw1); 176 ierr = MPI_Allreduce(w1,rw1,2*size,MPI_INT,PetscMaxSum_Op,comm);CHKERRQ(ierr); 177 bsz = rw1[rank]; 178 nrqr = rw1[size+rank]; 179 ierr = PetscFree(rw1);CHKERRQ(ierr); 180 } 181 182 /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */ 183 len = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int); 184 rbuf1 = (int**)PetscMalloc(len);CHKPTRQ(rbuf1); 185 rbuf1[0] = (int*)(rbuf1 + nrqr); 186 for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz; 187 188 /* Post the receives */ 189 r_waits1 = (MPI_Request*)PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1); 190 for (i=0; i<nrqr; ++i) { 191 ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr); 192 } 193 194 /* Allocate Memory for outgoing messages */ 195 len = 2*size*sizeof(int*) + 2*msz*sizeof(int)+ size*sizeof(int); 196 sbuf1 = (int **)PetscMalloc(len);CHKPTRQ(sbuf1); 197 ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 198 ierr = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr); 199 /* allocate memory for outgoing data + buf to receive the first reply */ 200 tmp = (int*)(ptr + size); 201 ctr = tmp + 2*msz; 202 203 { 204 int *iptr = tmp,ict = 0; 205 for (i=0; i<nrqs; i++) { 206 j = pa[i]; 207 iptr += ict; 208 sbuf1[j] = iptr; 209 ict = w1[j]; 210 } 211 } 212 213 /* Form the outgoing messages */ 214 /* Initialize the header space */ 215 for (i=0; i<nrqs; i++) { 216 j = pa[i]; 217 sbuf1[j][0] = 0; 218 ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr); 219 ptr[j] = sbuf1[j] + 2*w3[j] + 1; 220 } 221 222 /* Parse the isrow and copy data into outbuf */ 223 for (i=0; i<ismax; i++) { 224 ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 225 irow_i = irow[i]; 226 jmax = nrow[i]; 227 for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 228 row = irow_i[j]; 229 proc = rtable[row]; 230 if (proc != rank) { /* copy to the outgoing buf*/ 231 ctr[proc]++; 232 *ptr[proc] = row; 233 ptr[proc]++; 234 } 235 } 236 /* Update the headers for the current IS */ 237 for (j=0; j<size; j++) { /* Can Optimise this loop too */ 238 if ((ctr_j = ctr[j])) { 239 sbuf1_j = sbuf1[j]; 240 k = ++sbuf1_j[0]; 241 sbuf1_j[2*k] = ctr_j; 242 sbuf1_j[2*k-1] = i; 243 } 244 } 245 } 246 247 /* Now post the sends */ 248 s_waits1 = (MPI_Request*)PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1); 249 for (i=0; i<nrqs; ++i) { 250 j = pa[i]; 251 ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 252 } 253 254 /* Post recieves to capture the row_data from other procs */ 255 r_waits2 = (MPI_Request*)PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2); 256 rbuf2 = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar *));CHKPTRQ(rbuf2); 257 for (i=0; i<nrqs; i++) { 258 j = pa[i]; 259 count = (w1[j] - (2*sbuf1[j][0] + 1))*N; 260 rbuf2[i] = (Scalar *)PetscMalloc((count+1)*sizeof(Scalar));CHKPTRQ(rbuf2[i]); 261 ierr = MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 262 } 263 264 /* Receive messages(row_nos) and then, pack and send off the rowvalues 265 to the correct processors */ 266 267 s_waits2 = (MPI_Request*)PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2); 268 r_status1 = (MPI_Status*)PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(r_status1); 269 sbuf2 = (Scalar**)PetscMalloc((nrqr+1)*sizeof(Scalar*));CHKPTRQ(sbuf2); 270 271 { 272 Scalar *sbuf2_i,*v_start; 273 int s_proc; 274 for (i=0; i<nrqr; ++i) { 275 ierr = MPI_Waitany(nrqr,r_waits1,&index,r_status1+i);CHKERRQ(ierr); 276 s_proc = r_status1[i].MPI_SOURCE; /* send processor */ 277 rbuf1_i = rbuf1[index]; /* Actual message from s_proc */ 278 /* no of rows = end - start; since start is array index[], 0index, whel end 279 is length of the buffer - which is 1index */ 280 start = 2*rbuf1_i[0] + 1; 281 MPI_Get_count(r_status1+i,MPI_INT,&end); 282 /* allocate memory sufficinet to hold all the row values */ 283 sbuf2[index] = (Scalar *)PetscMalloc((end-start)*N*sizeof(Scalar));CHKPTRQ(sbuf2[index]); 284 sbuf2_i = sbuf2[index]; 285 /* Now pack the data */ 286 for (j=start; j<end; j++) { 287 row = rbuf1_i[j] - rstart; 288 v_start = a->v + row; 289 for (k=0; k<N; k++) { 290 sbuf2_i[0] = v_start[0]; 291 sbuf2_i++; v_start+=a->m; 292 } 293 } 294 /* Now send off the data */ 295 ierr = MPI_Isend(sbuf2[index],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);CHKERRQ(ierr); 296 } 297 } 298 /* End Send-Recv of IS + row_numbers */ 299 ierr = PetscFree(r_status1);CHKERRQ(ierr); 300 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 301 s_status1 = (MPI_Status*)PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status1); 302 ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr); 303 ierr = PetscFree(s_status1);CHKERRQ(ierr); 304 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 305 306 /* Create the submatrices */ 307 if (scall == MAT_REUSE_MATRIX) { 308 for (i=0; i<ismax; i++) { 309 mat = (Mat_SeqDense *)(submats[i]->data); 310 if ((mat->m != nrow[i]) || (mat->n != ncol[i])) { 311 SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size"); 312 } 313 ierr = PetscMemzero(mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 314 submats[i]->factor = C->factor; 315 } 316 } else { 317 for (i=0; i<ismax; i++) { 318 ierr = MatCreateSeqDense(PETSC_COMM_SELF,nrow[i],ncol[i],PETSC_NULL,submats+i);CHKERRQ(ierr); 319 } 320 } 321 322 /* Assemble the matrices */ 323 { 324 int col; 325 Scalar *imat_v,*mat_v,*imat_vi,*mat_vi; 326 327 for (i=0; i<ismax; i++) { 328 mat = (Mat_SeqDense*)submats[i]->data; 329 mat_v = a->v; 330 imat_v = mat->v; 331 irow_i = irow[i]; 332 m = nrow[i]; 333 for (j=0; j<m; j++) { 334 row = irow_i[j] ; 335 proc = rtable[row]; 336 if (proc == rank) { 337 row = row - rstart; 338 mat_vi = mat_v + row; 339 imat_vi = imat_v + j; 340 for (k=0; k<ncol[i]; k++) { 341 col = icol[i][k]; 342 imat_vi[k*m] = mat_vi[col*a->m]; 343 } 344 } 345 } 346 } 347 } 348 349 /* Create row map. This maps c->row to submat->row for each submat*/ 350 /* this is a very expensive operation wrt memory usage */ 351 len = (1+ismax)*sizeof(int*)+ ismax*c->M*sizeof(int); 352 rmap = (int **)PetscMalloc(len);CHKPTRQ(rmap); 353 rmap[0] = (int *)(rmap + ismax); 354 ierr = PetscMemzero(rmap[0],ismax*c->M*sizeof(int));CHKERRQ(ierr); 355 for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + c->M;} 356 for (i=0; i<ismax; i++) { 357 rmap_i = rmap[i]; 358 irow_i = irow[i]; 359 jmax = nrow[i]; 360 for (j=0; j<jmax; j++) { 361 rmap_i[irow_i[j]] = j; 362 } 363 } 364 365 /* Now Receive the row_values and assemble the rest of the matrix */ 366 367 r_status2 = (MPI_Status*)PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status2); 368 369 { 370 int is_max,tmp1,col,*sbuf1_i,is_sz; 371 Scalar *rbuf2_i,*imat_v,*imat_vi; 372 373 for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */ 374 ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr); 375 /* Now dig out the corresponding sbuf1, which contains the IS data_structure */ 376 sbuf1_i = sbuf1[pa[i]]; 377 is_max = sbuf1_i[0]; 378 ct1 = 2*is_max+1; 379 rbuf2_i = rbuf2[i]; 380 for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */ 381 is_no = sbuf1_i[2*j-1]; 382 is_sz = sbuf1_i[2*j]; 383 mat = (Mat_SeqDense*)submats[is_no]->data; 384 imat_v = mat->v; 385 rmap_i = rmap[is_no]; 386 m = nrow[is_no]; 387 for (k=0; k<is_sz; k++,rbuf2_i+=N) { /* For each row */ 388 row = sbuf1_i[ct1]; ct1++; 389 row = rmap_i[row]; 390 imat_vi = imat_v + row; 391 for (l=0; l<ncol[is_no]; l++) { /* For each col */ 392 col = icol[is_no][l]; 393 imat_vi[l*m] = rbuf2_i[col]; 394 } 395 } 396 } 397 } 398 } 399 /* End Send-Recv of row_values */ 400 ierr = PetscFree(r_status2);CHKERRQ(ierr); 401 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 402 s_status2 = (MPI_Status*)PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status1); 403 ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr); 404 ierr = PetscFree(s_status2);CHKERRQ(ierr); 405 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 406 407 /* Restore the indices */ 408 for (i=0; i<ismax; i++) { 409 ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 410 ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 411 } 412 413 /* Destroy allocated memory */ 414 ierr = PetscFree(irow);CHKERRQ(ierr); 415 ierr = PetscFree(w1);CHKERRQ(ierr); 416 ierr = PetscFree(pa);CHKERRQ(ierr); 417 418 419 for (i=0; i<nrqs; ++i) { 420 ierr = PetscFree(rbuf2[i]);CHKERRQ(ierr); 421 } 422 ierr = PetscFree(rbuf2);CHKERRQ(ierr); 423 ierr = PetscFree(sbuf1);CHKERRQ(ierr); 424 ierr = PetscFree(rbuf1);CHKERRQ(ierr); 425 426 for (i=0; i<nrqr; ++i) { 427 ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 428 } 429 430 ierr = PetscFree(sbuf2);CHKERRQ(ierr); 431 ierr = PetscFree(rmap);CHKERRQ(ierr); 432 433 ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 434 435 for (i=0; i<ismax; i++) { 436 ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 437 ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 438 } 439 440 PetscFunctionReturn(0); 441 } 442 443