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