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