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