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