1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: mpidense.c,v 1.98 1998/10/01 18:54:41 bsmith Exp bsmith $"; 3 #endif 4 5 /* 6 Basic functions for basic parallel dense matrices. 7 */ 8 9 #include "src/mat/impls/dense/mpi/mpidense.h" 10 #include "src/vec/vecimpl.h" 11 12 #undef __FUNC__ 13 #define __FUNC__ "MatSetValues_MPIDense" 14 int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 15 { 16 Mat_MPIDense *A = (Mat_MPIDense *) mat->data; 17 int ierr, i, j, rstart = A->rstart, rend = A->rend, row; 18 int roworiented = A->roworiented; 19 20 PetscFunctionBegin; 21 for ( i=0; i<m; i++ ) { 22 if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 23 if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 24 if (idxm[i] >= rstart && idxm[i] < rend) { 25 row = idxm[i] - rstart; 26 if (roworiented) { 27 ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv); CHKERRQ(ierr); 28 } else { 29 for ( j=0; j<n; j++ ) { 30 if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 31 if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 32 ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv); CHKERRQ(ierr); 33 } 34 } 35 } else { 36 if (roworiented) { 37 ierr = StashValues_Private(&A->stash,idxm[i],n,idxn,v+i*n,addv); CHKERRQ(ierr); 38 } else { /* must stash each seperately */ 39 row = idxm[i]; 40 for ( j=0; j<n; j++ ) { 41 ierr = StashValues_Private(&A->stash,row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 42 } 43 } 44 } 45 } 46 PetscFunctionReturn(0); 47 } 48 49 #undef __FUNC__ 50 #define __FUNC__ "MatGetValues_MPIDense" 51 int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 52 { 53 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 54 int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 55 56 PetscFunctionBegin; 57 for ( i=0; i<m; i++ ) { 58 if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 59 if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 60 if (idxm[i] >= rstart && idxm[i] < rend) { 61 row = idxm[i] - rstart; 62 for ( j=0; j<n; j++ ) { 63 if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 64 if (idxn[j] >= mdn->N) { 65 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 66 } 67 ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j); CHKERRQ(ierr); 68 } 69 } else { 70 SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 71 } 72 } 73 PetscFunctionReturn(0); 74 } 75 76 #undef __FUNC__ 77 #define __FUNC__ "MatGetArray_MPIDense" 78 int MatGetArray_MPIDense(Mat A,Scalar **array) 79 { 80 Mat_MPIDense *a = (Mat_MPIDense *) A->data; 81 int ierr; 82 83 PetscFunctionBegin; 84 ierr = MatGetArray(a->A,array); CHKERRQ(ierr); 85 PetscFunctionReturn(0); 86 } 87 88 #undef __FUNC__ 89 #define __FUNC__ "MatRestoreArray_MPIDense" 90 int MatRestoreArray_MPIDense(Mat A,Scalar **array) 91 { 92 PetscFunctionBegin; 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNC__ 97 #define __FUNC__ "MatAssemblyBegin_MPIDense" 98 int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 99 { 100 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 101 MPI_Comm comm = mat->comm; 102 int size = mdn->size, *owners = mdn->rowners, rank = mdn->rank; 103 int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 104 int tag = mat->tag, *owner,*starts,count,ierr; 105 InsertMode addv; 106 MPI_Request *send_waits,*recv_waits; 107 Scalar *rvalues,*svalues; 108 109 PetscFunctionBegin; 110 /* make sure all processors are either in INSERTMODE or ADDMODE */ 111 ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 112 if (addv == (ADD_VALUES|INSERT_VALUES)) { 113 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs"); 114 } 115 mat->insertmode = addv; /* in case this processor had no cache */ 116 117 /* first count number of contributors to each processor */ 118 nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 119 PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 120 owner = (int *) PetscMalloc( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 121 for ( i=0; i<mdn->stash.n; i++ ) { 122 idx = mdn->stash.idx[i]; 123 for ( j=0; j<size; j++ ) { 124 if (idx >= owners[j] && idx < owners[j+1]) { 125 nprocs[j]++; procs[j] = 1; owner[i] = j; break; 126 } 127 } 128 } 129 nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 130 131 /* inform other processors of number of messages and max length*/ 132 work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 133 ierr = MPI_Allreduce(procs,work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 134 nreceives = work[rank]; 135 if (nreceives > size) SETERRQ(PETSC_ERR_PLIB,0,"Internal PETSc error"); 136 ierr = MPI_Allreduce(nprocs,work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 137 nmax = work[rank]; 138 PetscFree(work); 139 140 /* post receives: 141 1) each message will consist of ordered pairs 142 (global index,value) we store the global index as a double 143 to simplify the message passing. 144 2) since we don't know how long each individual message is we 145 allocate the largest needed buffer for each receive. Potentially 146 this is a lot of wasted space. 147 148 This could be done better. 149 */ 150 rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));CHKPTRQ(rvalues); 151 recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 152 for ( i=0; i<nreceives; i++ ) { 153 ierr = MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 154 } 155 156 /* do sends: 157 1) starts[i] gives the starting index in svalues for stuff going to 158 the ith processor 159 */ 160 svalues = (Scalar *) PetscMalloc( 3*(mdn->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues); 161 send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 162 starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 163 starts[0] = 0; 164 for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 165 for ( i=0; i<mdn->stash.n; i++ ) { 166 svalues[3*starts[owner[i]]] = (Scalar) mdn->stash.idx[i]; 167 svalues[3*starts[owner[i]]+1] = (Scalar) mdn->stash.idy[i]; 168 svalues[3*(starts[owner[i]]++)+2] = mdn->stash.array[i]; 169 } 170 PetscFree(owner); 171 starts[0] = 0; 172 for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 173 count = 0; 174 for ( i=0; i<size; i++ ) { 175 if (procs[i]) { 176 ierr = MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 177 } 178 } 179 PetscFree(starts); PetscFree(nprocs); 180 181 /* Free cache space */ 182 PLogInfo(mat,"MatAssemblyBegin_MPIDense:Number of off-processor values %d\n",mdn->stash.n); 183 ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr); 184 185 mdn->svalues = svalues; mdn->rvalues = rvalues; 186 mdn->nsends = nsends; mdn->nrecvs = nreceives; 187 mdn->send_waits = send_waits; mdn->recv_waits = recv_waits; 188 mdn->rmax = nmax; 189 190 PetscFunctionReturn(0); 191 } 192 extern int MatSetUpMultiply_MPIDense(Mat); 193 194 #undef __FUNC__ 195 #define __FUNC__ "MatAssemblyEnd_MPIDense" 196 int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 197 { 198 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 199 MPI_Status *send_status,recv_status; 200 int imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col; 201 Scalar *values,val; 202 InsertMode addv = mat->insertmode; 203 204 PetscFunctionBegin; 205 /* wait on receives */ 206 while (count) { 207 ierr = MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 208 /* unpack receives into our local space */ 209 values = mdn->rvalues + 3*imdex*mdn->rmax; 210 ierr = MPI_Get_count(&recv_status,MPIU_SCALAR,&n);CHKERRQ(ierr); 211 n = n/3; 212 for ( i=0; i<n; i++ ) { 213 row = (int) PetscReal(values[3*i]) - mdn->rstart; 214 col = (int) PetscReal(values[3*i+1]); 215 val = values[3*i+2]; 216 if (col >= 0 && col < mdn->N) { 217 MatSetValues(mdn->A,1,&row,1,&col,&val,addv); 218 } 219 else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid column");} 220 } 221 count--; 222 } 223 PetscFree(mdn->recv_waits); PetscFree(mdn->rvalues); 224 225 /* wait on sends */ 226 if (mdn->nsends) { 227 send_status = (MPI_Status *) PetscMalloc(mdn->nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 228 ierr = MPI_Waitall(mdn->nsends,mdn->send_waits,send_status);CHKERRQ(ierr); 229 PetscFree(send_status); 230 } 231 PetscFree(mdn->send_waits); PetscFree(mdn->svalues); 232 233 ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr); 234 ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr); 235 236 if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 237 ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr); 238 } 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNC__ 243 #define __FUNC__ "MatZeroEntries_MPIDense" 244 int MatZeroEntries_MPIDense(Mat A) 245 { 246 int ierr; 247 Mat_MPIDense *l = (Mat_MPIDense *) A->data; 248 249 PetscFunctionBegin; 250 ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNC__ 255 #define __FUNC__ "MatGetBlockSize_MPIDense" 256 int MatGetBlockSize_MPIDense(Mat A,int *bs) 257 { 258 PetscFunctionBegin; 259 *bs = 1; 260 PetscFunctionReturn(0); 261 } 262 263 /* the code does not do the diagonal entries correctly unless the 264 matrix is square and the column and row owerships are identical. 265 This is a BUG. The only way to fix it seems to be to access 266 mdn->A and mdn->B directly and not through the MatZeroRows() 267 routine. 268 */ 269 #undef __FUNC__ 270 #define __FUNC__ "MatZeroRows_MPIDense" 271 int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 272 { 273 Mat_MPIDense *l = (Mat_MPIDense *) A->data; 274 int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 275 int *procs,*nprocs,j,found,idx,nsends,*work; 276 int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 277 int *rvalues,tag = A->tag,count,base,slen,n,*source; 278 int *lens,imdex,*lrows,*values; 279 MPI_Comm comm = A->comm; 280 MPI_Request *send_waits,*recv_waits; 281 MPI_Status recv_status,*send_status; 282 IS istmp; 283 284 PetscFunctionBegin; 285 ierr = ISGetSize(is,&N); CHKERRQ(ierr); 286 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 287 288 /* first count number of contributors to each processor */ 289 nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 290 PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 291 owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 292 for ( i=0; i<N; i++ ) { 293 idx = rows[i]; 294 found = 0; 295 for ( j=0; j<size; j++ ) { 296 if (idx >= owners[j] && idx < owners[j+1]) { 297 nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 298 } 299 } 300 if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 301 } 302 nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 303 304 /* inform other processors of number of messages and max length*/ 305 work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 306 ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 307 nrecvs = work[rank]; 308 ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 309 nmax = work[rank]; 310 PetscFree(work); 311 312 /* post receives: */ 313 rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 314 recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 315 for ( i=0; i<nrecvs; i++ ) { 316 ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 317 } 318 319 /* do sends: 320 1) starts[i] gives the starting index in svalues for stuff going to 321 the ith processor 322 */ 323 svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 324 send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 325 starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 326 starts[0] = 0; 327 for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 328 for ( i=0; i<N; i++ ) { 329 svalues[starts[owner[i]]++] = rows[i]; 330 } 331 ISRestoreIndices(is,&rows); 332 333 starts[0] = 0; 334 for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 335 count = 0; 336 for ( i=0; i<size; i++ ) { 337 if (procs[i]) { 338 ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 339 } 340 } 341 PetscFree(starts); 342 343 base = owners[rank]; 344 345 /* wait on receives */ 346 lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 347 source = lens + nrecvs; 348 count = nrecvs; slen = 0; 349 while (count) { 350 ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 351 /* unpack receives into our local space */ 352 ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 353 source[imdex] = recv_status.MPI_SOURCE; 354 lens[imdex] = n; 355 slen += n; 356 count--; 357 } 358 PetscFree(recv_waits); 359 360 /* move the data into the send scatter */ 361 lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 362 count = 0; 363 for ( i=0; i<nrecvs; i++ ) { 364 values = rvalues + i*nmax; 365 for ( j=0; j<lens[i]; j++ ) { 366 lrows[count++] = values[j] - base; 367 } 368 } 369 PetscFree(rvalues); PetscFree(lens); 370 PetscFree(owner); PetscFree(nprocs); 371 372 /* actually zap the local rows */ 373 ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 374 PLogObjectParent(A,istmp); 375 PetscFree(lrows); 376 ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 377 ierr = ISDestroy(istmp); CHKERRQ(ierr); 378 379 /* wait on sends */ 380 if (nsends) { 381 send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 382 ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 383 PetscFree(send_status); 384 } 385 PetscFree(send_waits); PetscFree(svalues); 386 387 PetscFunctionReturn(0); 388 } 389 390 #undef __FUNC__ 391 #define __FUNC__ "MatMult_MPIDense" 392 int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 393 { 394 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 395 int ierr; 396 397 PetscFunctionBegin; 398 ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 399 ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 400 ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy); CHKERRQ(ierr); 401 PetscFunctionReturn(0); 402 } 403 404 #undef __FUNC__ 405 #define __FUNC__ "MatMultAdd_MPIDense" 406 int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 407 { 408 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 409 int ierr; 410 411 PetscFunctionBegin; 412 ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 413 ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 414 ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr); 415 PetscFunctionReturn(0); 416 } 417 418 #undef __FUNC__ 419 #define __FUNC__ "MatMultTrans_MPIDense" 420 int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 421 { 422 Mat_MPIDense *a = (Mat_MPIDense *) A->data; 423 int ierr; 424 Scalar zero = 0.0; 425 426 PetscFunctionBegin; 427 ierr = VecSet(&zero,yy); CHKERRQ(ierr); 428 ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr); 429 ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 430 ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNC__ 435 #define __FUNC__ "MatMultTransAdd_MPIDense" 436 int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 437 { 438 Mat_MPIDense *a = (Mat_MPIDense *) A->data; 439 int ierr; 440 441 PetscFunctionBegin; 442 ierr = VecCopy(yy,zz); CHKERRQ(ierr); 443 ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr); 444 ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 445 ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 446 PetscFunctionReturn(0); 447 } 448 449 #undef __FUNC__ 450 #define __FUNC__ "MatGetDiagonal_MPIDense" 451 int MatGetDiagonal_MPIDense(Mat A,Vec v) 452 { 453 Mat_MPIDense *a = (Mat_MPIDense *) A->data; 454 Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 455 int ierr, len, i, n, m = a->m, radd; 456 Scalar *x, zero = 0.0; 457 458 PetscFunctionBegin; 459 VecSet(&zero,v); 460 ierr = VecGetArray(v,&x); CHKERRQ(ierr); 461 ierr = VecGetSize(v,&n); CHKERRQ(ierr); 462 if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 463 len = PetscMin(aloc->m,aloc->n); 464 radd = a->rstart*m; 465 for ( i=0; i<len; i++ ) { 466 x[i] = aloc->v[radd + i*m + i]; 467 } 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNC__ 472 #define __FUNC__ "MatDestroy_MPIDense" 473 int MatDestroy_MPIDense(Mat mat) 474 { 475 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 476 int ierr; 477 478 PetscFunctionBegin; 479 if (--mat->refct > 0) PetscFunctionReturn(0); 480 481 if (mat->mapping) { 482 ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 483 } 484 if (mat->bmapping) { 485 ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr); 486 } 487 #if defined(USE_PETSC_LOG) 488 PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N); 489 #endif 490 PetscFree(mdn->rowners); 491 ierr = MatDestroy(mdn->A); CHKERRQ(ierr); 492 if (mdn->lvec) VecDestroy(mdn->lvec); 493 if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 494 if (mdn->factor) { 495 if (mdn->factor->temp) PetscFree(mdn->factor->temp); 496 if (mdn->factor->tag) PetscFree(mdn->factor->tag); 497 if (mdn->factor->pivots) PetscFree(mdn->factor->pivots); 498 PetscFree(mdn->factor); 499 } 500 PetscFree(mdn); 501 if (mat->rmap) { 502 ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 503 } 504 if (mat->cmap) { 505 ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 506 } 507 PLogObjectDestroy(mat); 508 PetscHeaderDestroy(mat); 509 PetscFunctionReturn(0); 510 } 511 512 #undef __FUNC__ 513 #define __FUNC__ "MatView_MPIDense_Binary" 514 static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 515 { 516 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 517 int ierr; 518 519 PetscFunctionBegin; 520 if (mdn->size == 1) { 521 ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 522 } 523 else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 524 PetscFunctionReturn(0); 525 } 526 527 #undef __FUNC__ 528 #define __FUNC__ "MatView_MPIDense_ASCII" 529 static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 530 { 531 Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 532 int ierr, format, size = mdn->size, rank = mdn->rank; 533 FILE *fd; 534 ViewerType vtype; 535 536 PetscFunctionBegin; 537 ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 538 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 539 ierr = ViewerGetFormat(viewer,&format); 540 if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 541 int rank; 542 MatInfo info; 543 MPI_Comm_rank(mat->comm,&rank); 544 ierr = MatGetInfo(mat,MAT_LOCAL,&info); 545 PetscSequentialPhaseBegin(mat->comm,1); 546 fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m, 547 (int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 548 fflush(fd); 549 PetscSequentialPhaseEnd(mat->comm,1); 550 ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr); 551 PetscFunctionReturn(0); 552 } 553 else if (format == VIEWER_FORMAT_ASCII_INFO) { 554 PetscFunctionReturn(0); 555 } 556 557 if (size == 1) { 558 ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 559 } else { 560 /* assemble the entire matrix onto first processor. */ 561 Mat A; 562 int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 563 Scalar *vals; 564 Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 565 566 if (!rank) { 567 ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr); 568 } else { 569 ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr); 570 } 571 PLogObjectParent(mat,A); 572 573 /* Copy the matrix ... This isn't the most efficient means, 574 but it's quick for now */ 575 row = mdn->rstart; m = Amdn->m; 576 for ( i=0; i<m; i++ ) { 577 ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 578 ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 579 ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 580 row++; 581 } 582 583 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 584 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 585 if (!rank) { 586 ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr); 587 } 588 ierr = MatDestroy(A); CHKERRQ(ierr); 589 } 590 PetscFunctionReturn(0); 591 } 592 593 #undef __FUNC__ 594 #define __FUNC__ "MatView_MPIDense" 595 int MatView_MPIDense(Mat mat,Viewer viewer) 596 { 597 int ierr; 598 ViewerType vtype; 599 600 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 601 if (!PetscStrcmp(vtype,ASCII_VIEWER)) { 602 ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 603 } else if (!PetscStrcmp(vtype,BINARY_VIEWER)) { 604 ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 605 } else { 606 SETERRQ(1,1,"Viewer type not supported by PETSc object"); 607 } 608 PetscFunctionReturn(0); 609 } 610 611 #undef __FUNC__ 612 #define __FUNC__ "MatGetInfo_MPIDense" 613 int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 614 { 615 Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 616 Mat mdn = mat->A; 617 int ierr; 618 double isend[5], irecv[5]; 619 620 PetscFunctionBegin; 621 info->rows_global = (double)mat->M; 622 info->columns_global = (double)mat->N; 623 info->rows_local = (double)mat->m; 624 info->columns_local = (double)mat->N; 625 info->block_size = 1.0; 626 ierr = MatGetInfo(mdn,MAT_LOCAL,info); CHKERRQ(ierr); 627 isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 628 isend[3] = info->memory; isend[4] = info->mallocs; 629 if (flag == MAT_LOCAL) { 630 info->nz_used = isend[0]; 631 info->nz_allocated = isend[1]; 632 info->nz_unneeded = isend[2]; 633 info->memory = isend[3]; 634 info->mallocs = isend[4]; 635 } else if (flag == MAT_GLOBAL_MAX) { 636 ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 637 info->nz_used = irecv[0]; 638 info->nz_allocated = irecv[1]; 639 info->nz_unneeded = irecv[2]; 640 info->memory = irecv[3]; 641 info->mallocs = irecv[4]; 642 } else if (flag == MAT_GLOBAL_SUM) { 643 ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 644 info->nz_used = irecv[0]; 645 info->nz_allocated = irecv[1]; 646 info->nz_unneeded = irecv[2]; 647 info->memory = irecv[3]; 648 info->mallocs = irecv[4]; 649 } 650 info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 651 info->fill_ratio_needed = 0; 652 info->factor_mallocs = 0; 653 PetscFunctionReturn(0); 654 } 655 656 /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*); 657 extern int MatLUFactorNumeric_MPIDense(Mat,Mat*); 658 extern int MatLUFactor_MPIDense(Mat,IS,IS,double); 659 extern int MatSolve_MPIDense(Mat,Vec,Vec); 660 extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec); 661 extern int MatSolveTrans_MPIDense(Mat,Vec,Vec); 662 extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */ 663 664 #undef __FUNC__ 665 #define __FUNC__ "MatSetOption_MPIDense" 666 int MatSetOption_MPIDense(Mat A,MatOption op) 667 { 668 Mat_MPIDense *a = (Mat_MPIDense *) A->data; 669 670 PetscFunctionBegin; 671 if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 672 op == MAT_YES_NEW_NONZERO_LOCATIONS || 673 op == MAT_NEW_NONZERO_LOCATION_ERROR || 674 op == MAT_NEW_NONZERO_ALLOCATION_ERROR || 675 op == MAT_COLUMNS_SORTED || 676 op == MAT_COLUMNS_UNSORTED) { 677 MatSetOption(a->A,op); 678 } else if (op == MAT_ROW_ORIENTED) { 679 a->roworiented = 1; 680 MatSetOption(a->A,op); 681 } else if (op == MAT_ROWS_SORTED || 682 op == MAT_ROWS_UNSORTED || 683 op == MAT_SYMMETRIC || 684 op == MAT_STRUCTURALLY_SYMMETRIC || 685 op == MAT_YES_NEW_DIAGONALS || 686 op == MAT_USE_HASH_TABLE) { 687 PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 688 } else if (op == MAT_COLUMN_ORIENTED) { 689 a->roworiented = 0; MatSetOption(a->A,op); 690 } else if (op == MAT_NO_NEW_DIAGONALS) { 691 SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 692 } else { 693 SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 694 } 695 PetscFunctionReturn(0); 696 } 697 698 #undef __FUNC__ 699 #define __FUNC__ "MatGetSize_MPIDense" 700 int MatGetSize_MPIDense(Mat A,int *m,int *n) 701 { 702 Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 703 704 PetscFunctionBegin; 705 *m = mat->M; *n = mat->N; 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNC__ 710 #define __FUNC__ "MatGetLocalSize_MPIDense" 711 int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 712 { 713 Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 714 715 PetscFunctionBegin; 716 *m = mat->m; *n = mat->N; 717 PetscFunctionReturn(0); 718 } 719 720 #undef __FUNC__ 721 #define __FUNC__ "MatGetOwnershipRange_MPIDense" 722 int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 723 { 724 Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 725 726 PetscFunctionBegin; 727 *m = mat->rstart; *n = mat->rend; 728 PetscFunctionReturn(0); 729 } 730 731 #undef __FUNC__ 732 #define __FUNC__ "MatGetRow_MPIDense" 733 int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 734 { 735 Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 736 int lrow, rstart = mat->rstart, rend = mat->rend,ierr; 737 738 PetscFunctionBegin; 739 if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows") 740 lrow = row - rstart; 741 ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 742 PetscFunctionReturn(0); 743 } 744 745 #undef __FUNC__ 746 #define __FUNC__ "MatRestoreRow_MPIDense" 747 int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 748 { 749 PetscFunctionBegin; 750 if (idx) PetscFree(*idx); 751 if (v) PetscFree(*v); 752 PetscFunctionReturn(0); 753 } 754 755 #undef __FUNC__ 756 #define __FUNC__ "MatNorm_MPIDense" 757 int MatNorm_MPIDense(Mat A,NormType type,double *norm) 758 { 759 Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 760 Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 761 int ierr, i, j; 762 double sum = 0.0; 763 Scalar *v = mat->v; 764 765 PetscFunctionBegin; 766 if (mdn->size == 1) { 767 ierr = MatNorm(mdn->A,type,norm); CHKERRQ(ierr); 768 } else { 769 if (type == NORM_FROBENIUS) { 770 for (i=0; i<mat->n*mat->m; i++ ) { 771 #if defined(USE_PETSC_COMPLEX) 772 sum += PetscReal(PetscConj(*v)*(*v)); v++; 773 #else 774 sum += (*v)*(*v); v++; 775 #endif 776 } 777 ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 778 *norm = sqrt(*norm); 779 PLogFlops(2*mat->n*mat->m); 780 } else if (type == NORM_1) { 781 double *tmp, *tmp2; 782 tmp = (double *) PetscMalloc( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp); 783 tmp2 = tmp + mdn->N; 784 PetscMemzero(tmp,2*mdn->N*sizeof(double)); 785 *norm = 0.0; 786 v = mat->v; 787 for ( j=0; j<mat->n; j++ ) { 788 for ( i=0; i<mat->m; i++ ) { 789 tmp[j] += PetscAbsScalar(*v); v++; 790 } 791 } 792 ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 793 for ( j=0; j<mdn->N; j++ ) { 794 if (tmp2[j] > *norm) *norm = tmp2[j]; 795 } 796 PetscFree(tmp); 797 PLogFlops(mat->n*mat->m); 798 } else if (type == NORM_INFINITY) { /* max row norm */ 799 double ntemp; 800 ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr); 801 ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 802 } else { 803 SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 804 } 805 } 806 PetscFunctionReturn(0); 807 } 808 809 #undef __FUNC__ 810 #define __FUNC__ "MatTranspose_MPIDense" 811 int MatTranspose_MPIDense(Mat A,Mat *matout) 812 { 813 Mat_MPIDense *a = (Mat_MPIDense *) A->data; 814 Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 815 Mat B; 816 int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 817 int j, i, ierr; 818 Scalar *v; 819 820 PetscFunctionBegin; 821 if (matout == PETSC_NULL && M != N) { 822 SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 823 } 824 ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 825 826 m = Aloc->m; n = Aloc->n; v = Aloc->v; 827 rwork = (int *) PetscMalloc(n*sizeof(int)); CHKPTRQ(rwork); 828 for ( j=0; j<n; j++ ) { 829 for (i=0; i<m; i++) rwork[i] = rstart + i; 830 ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr); 831 v += m; 832 } 833 PetscFree(rwork); 834 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 835 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 836 if (matout != PETSC_NULL) { 837 *matout = B; 838 } else { 839 PetscOps *Abops; 840 MatOps Aops; 841 842 /* This isn't really an in-place transpose, but free data struct from a */ 843 PetscFree(a->rowners); 844 ierr = MatDestroy(a->A); CHKERRQ(ierr); 845 if (a->lvec) VecDestroy(a->lvec); 846 if (a->Mvctx) VecScatterDestroy(a->Mvctx); 847 PetscFree(a); 848 849 /* 850 This is horrible, horrible code. We need to keep the 851 A pointers for the bops and ops but copy everything 852 else from C. 853 */ 854 Abops = A->bops; 855 Aops = A->ops; 856 PetscMemcpy(A,B,sizeof(struct _p_Mat)); 857 A->bops = Abops; 858 A->ops = Aops; 859 860 PetscHeaderDestroy(B); 861 } 862 PetscFunctionReturn(0); 863 } 864 865 #include "pinclude/blaslapack.h" 866 #undef __FUNC__ 867 #define __FUNC__ "MatScale_MPIDense" 868 int MatScale_MPIDense(Scalar *alpha,Mat inA) 869 { 870 Mat_MPIDense *A = (Mat_MPIDense *) inA->data; 871 Mat_SeqDense *a = (Mat_SeqDense *) A->A->data; 872 int one = 1, nz; 873 874 PetscFunctionBegin; 875 nz = a->m*a->n; 876 BLscal_( &nz, alpha, a->v, &one ); 877 PLogFlops(nz); 878 PetscFunctionReturn(0); 879 } 880 881 static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 882 extern int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatGetSubMatrixCall,Mat **); 883 884 /* -------------------------------------------------------------------*/ 885 static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 886 MatGetRow_MPIDense, 887 MatRestoreRow_MPIDense, 888 MatMult_MPIDense, 889 MatMultAdd_MPIDense, 890 MatMultTrans_MPIDense, 891 MatMultTransAdd_MPIDense, 892 0, 893 0, 894 0, 895 0, 896 0, 897 0, 898 0, 899 MatTranspose_MPIDense, 900 MatGetInfo_MPIDense,0, 901 MatGetDiagonal_MPIDense, 902 0, 903 MatNorm_MPIDense, 904 MatAssemblyBegin_MPIDense, 905 MatAssemblyEnd_MPIDense, 906 0, 907 MatSetOption_MPIDense, 908 MatZeroEntries_MPIDense, 909 MatZeroRows_MPIDense, 910 0, 911 0, 912 0, 913 0, 914 MatGetSize_MPIDense, 915 MatGetLocalSize_MPIDense, 916 MatGetOwnershipRange_MPIDense, 917 0, 918 0, 919 MatGetArray_MPIDense, 920 MatRestoreArray_MPIDense, 921 MatDuplicate_MPIDense, 922 0, 923 0, 924 0, 925 0, 926 0, 927 MatGetSubMatrices_MPIDense, 928 0, 929 MatGetValues_MPIDense, 930 0, 931 0, 932 MatScale_MPIDense, 933 0, 934 0, 935 0, 936 MatGetBlockSize_MPIDense, 937 0, 938 0, 939 0, 940 0, 941 0, 942 0, 943 0, 944 0, 945 0, 946 0, 947 0, 948 0, 949 MatGetMaps_Petsc}; 950 951 #undef __FUNC__ 952 #define __FUNC__ "MatCreateMPIDense" 953 /*@C 954 MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 955 956 Collective on MPI_Comm 957 958 Input Parameters: 959 + comm - MPI communicator 960 . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 961 . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 962 . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 963 . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 964 - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 965 to control all matrix memory allocation. 966 967 Output Parameter: 968 . A - the matrix 969 970 Notes: 971 The dense format is fully compatible with standard Fortran 77 972 storage by columns. 973 974 The data input variable is intended primarily for Fortran programmers 975 who wish to allocate their own matrix memory space. Most users should 976 set data=PETSC_NULL. 977 978 The user MUST specify either the local or global matrix dimensions 979 (possibly both). 980 981 Currently, the only parallel dense matrix decomposition is by rows, 982 so that n=N and each submatrix owns all of the global columns. 983 984 .keywords: matrix, dense, parallel 985 986 .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 987 @*/ 988 int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 989 { 990 Mat mat; 991 Mat_MPIDense *a; 992 int ierr, i,flg; 993 994 PetscFunctionBegin; 995 /* Note: For now, when data is specified above, this assumes the user correctly 996 allocates the local dense storage space. We should add error checking. */ 997 998 *A = 0; 999 PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,comm,MatDestroy,MatView); 1000 PLogObjectCreate(mat); 1001 mat->data = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a); 1002 PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps)); 1003 mat->ops->destroy = MatDestroy_MPIDense; 1004 mat->ops->view = MatView_MPIDense; 1005 mat->factor = 0; 1006 mat->mapping = 0; 1007 1008 a->factor = 0; 1009 mat->insertmode = NOT_SET_VALUES; 1010 MPI_Comm_rank(comm,&a->rank); 1011 MPI_Comm_size(comm,&a->size); 1012 1013 if (M == PETSC_DECIDE) {ierr = MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);} 1014 if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);} 1015 1016 /* 1017 The computation of n is wrong below, n should represent the number of local 1018 rows in the right (column vector) 1019 */ 1020 1021 /* each row stores all columns */ 1022 if (N == PETSC_DECIDE) N = n; 1023 if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 1024 /* if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */ 1025 a->N = mat->N = N; 1026 a->M = mat->M = M; 1027 a->m = mat->m = m; 1028 a->n = mat->n = n; 1029 1030 /* the information in the maps duplicates the information computed below, eventually 1031 we should remove the duplicate information that is not contained in the maps */ 1032 ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr); 1033 ierr = MapCreateMPI(comm,n,N,&mat->cmap);CHKERRQ(ierr); 1034 1035 /* build local table of row and column ownerships */ 1036 a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 1037 a->cowners = a->rowners + a->size + 1; 1038 PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1039 ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1040 a->rowners[0] = 0; 1041 for ( i=2; i<=a->size; i++ ) { 1042 a->rowners[i] += a->rowners[i-1]; 1043 } 1044 a->rstart = a->rowners[a->rank]; 1045 a->rend = a->rowners[a->rank+1]; 1046 ierr = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1047 a->cowners[0] = 0; 1048 for ( i=2; i<=a->size; i++ ) { 1049 a->cowners[i] += a->cowners[i-1]; 1050 } 1051 1052 ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A); CHKERRQ(ierr); 1053 PLogObjectParent(mat,a->A); 1054 1055 /* build cache for off array entries formed */ 1056 ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr); 1057 1058 /* stuff used for matrix vector multiply */ 1059 a->lvec = 0; 1060 a->Mvctx = 0; 1061 a->roworiented = 1; 1062 1063 *A = mat; 1064 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 1065 if (flg) { 1066 ierr = MatPrintHelp(mat); CHKERRQ(ierr); 1067 } 1068 PetscFunctionReturn(0); 1069 } 1070 1071 #undef __FUNC__ 1072 #define __FUNC__ "MatDuplicate_MPIDense" 1073 static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 1074 { 1075 Mat mat; 1076 Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 1077 int ierr; 1078 FactorCtx *factor; 1079 1080 PetscFunctionBegin; 1081 *newmat = 0; 1082 PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,A->comm,MatDestroy,MatView); 1083 PLogObjectCreate(mat); 1084 mat->data = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a); 1085 PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps)); 1086 mat->ops->destroy = MatDestroy_MPIDense; 1087 mat->ops->view = MatView_MPIDense; 1088 mat->factor = A->factor; 1089 mat->assembled = PETSC_TRUE; 1090 1091 a->m = mat->m = oldmat->m; 1092 a->n = mat->n = oldmat->n; 1093 a->M = mat->M = oldmat->M; 1094 a->N = mat->N = oldmat->N; 1095 if (oldmat->factor) { 1096 a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx)); CHKPTRQ(factor); 1097 /* copy factor contents ... add this code! */ 1098 } else a->factor = 0; 1099 1100 a->rstart = oldmat->rstart; 1101 a->rend = oldmat->rend; 1102 a->size = oldmat->size; 1103 a->rank = oldmat->rank; 1104 mat->insertmode = NOT_SET_VALUES; 1105 1106 a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners); 1107 PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1108 PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int)); 1109 ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 1110 1111 ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 1112 PLogObjectParent(mat,a->lvec); 1113 ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 1114 PLogObjectParent(mat,a->Mvctx); 1115 ierr = MatDuplicate(oldmat->A,cpvalues,&a->A); CHKERRQ(ierr); 1116 PLogObjectParent(mat,a->A); 1117 *newmat = mat; 1118 PetscFunctionReturn(0); 1119 } 1120 1121 #include "sys.h" 1122 1123 #undef __FUNC__ 1124 #define __FUNC__ "MatLoad_MPIDense_DenseInFile" 1125 int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat) 1126 { 1127 int *rowners, i,size,rank,m,ierr,nz,j; 1128 Scalar *array,*vals,*vals_ptr; 1129 MPI_Status status; 1130 1131 PetscFunctionBegin; 1132 MPI_Comm_rank(comm,&rank); 1133 MPI_Comm_size(comm,&size); 1134 1135 /* determine ownership of all rows */ 1136 m = M/size + ((M % size) > rank); 1137 rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1138 ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1139 rowners[0] = 0; 1140 for ( i=2; i<=size; i++ ) { 1141 rowners[i] += rowners[i-1]; 1142 } 1143 1144 ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 1145 ierr = MatGetArray(*newmat,&array); CHKERRQ(ierr); 1146 1147 if (!rank) { 1148 vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals); 1149 1150 /* read in my part of the matrix numerical values */ 1151 ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR); CHKERRQ(ierr); 1152 1153 /* insert into matrix-by row (this is why cannot directly read into array */ 1154 vals_ptr = vals; 1155 for ( i=0; i<m; i++ ) { 1156 for ( j=0; j<N; j++ ) { 1157 array[i + j*m] = *vals_ptr++; 1158 } 1159 } 1160 1161 /* read in other processors and ship out */ 1162 for ( i=1; i<size; i++ ) { 1163 nz = (rowners[i+1] - rowners[i])*N; 1164 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1165 ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 1166 } 1167 } else { 1168 /* receive numeric values */ 1169 vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals); 1170 1171 /* receive message of values*/ 1172 ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 1173 1174 /* insert into matrix-by row (this is why cannot directly read into array */ 1175 vals_ptr = vals; 1176 for ( i=0; i<m; i++ ) { 1177 for ( j=0; j<N; j++ ) { 1178 array[i + j*m] = *vals_ptr++; 1179 } 1180 } 1181 } 1182 PetscFree(rowners); 1183 PetscFree(vals); 1184 ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1185 ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1186 PetscFunctionReturn(0); 1187 } 1188 1189 1190 #undef __FUNC__ 1191 #define __FUNC__ "MatLoad_MPIDense" 1192 int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 1193 { 1194 Mat A; 1195 Scalar *vals,*svals; 1196 MPI_Comm comm = ((PetscObject)viewer)->comm; 1197 MPI_Status status; 1198 int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1199 int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 1200 int tag = ((PetscObject)viewer)->tag; 1201 int i, nz, ierr, j,rstart, rend, fd; 1202 1203 PetscFunctionBegin; 1204 MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 1205 if (!rank) { 1206 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 1207 ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr); 1208 if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 1209 } 1210 1211 ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 1212 M = header[1]; N = header[2]; nz = header[3]; 1213 1214 /* 1215 Handle case where matrix is stored on disk as a dense matrix 1216 */ 1217 if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1218 ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 1219 PetscFunctionReturn(0); 1220 } 1221 1222 /* determine ownership of all rows */ 1223 m = M/size + ((M % size) > rank); 1224 rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1225 ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1226 rowners[0] = 0; 1227 for ( i=2; i<=size; i++ ) { 1228 rowners[i] += rowners[i-1]; 1229 } 1230 rstart = rowners[rank]; 1231 rend = rowners[rank+1]; 1232 1233 /* distribute row lengths to all processors */ 1234 ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 1235 offlens = ourlens + (rend-rstart); 1236 if (!rank) { 1237 rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1238 ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 1239 sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 1240 for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1241 ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1242 PetscFree(sndcounts); 1243 } else { 1244 ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 1245 } 1246 1247 if (!rank) { 1248 /* calculate the number of nonzeros on each processor */ 1249 procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1250 PetscMemzero(procsnz,size*sizeof(int)); 1251 for ( i=0; i<size; i++ ) { 1252 for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 1253 procsnz[i] += rowlengths[j]; 1254 } 1255 } 1256 PetscFree(rowlengths); 1257 1258 /* determine max buffer needed and allocate it */ 1259 maxnz = 0; 1260 for ( i=0; i<size; i++ ) { 1261 maxnz = PetscMax(maxnz,procsnz[i]); 1262 } 1263 cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 1264 1265 /* read in my part of the matrix column indices */ 1266 nz = procsnz[0]; 1267 mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1268 ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr); 1269 1270 /* read in every one elses and ship off */ 1271 for ( i=1; i<size; i++ ) { 1272 nz = procsnz[i]; 1273 ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr); 1274 ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 1275 } 1276 PetscFree(cols); 1277 } else { 1278 /* determine buffer space needed for message */ 1279 nz = 0; 1280 for ( i=0; i<m; i++ ) { 1281 nz += ourlens[i]; 1282 } 1283 mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 1284 1285 /* receive message of column indices*/ 1286 ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1287 ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1288 if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 1289 } 1290 1291 /* loop over local rows, determining number of off diagonal entries */ 1292 PetscMemzero(offlens,m*sizeof(int)); 1293 jj = 0; 1294 for ( i=0; i<m; i++ ) { 1295 for ( j=0; j<ourlens[i]; j++ ) { 1296 if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1297 jj++; 1298 } 1299 } 1300 1301 /* create our matrix */ 1302 for ( i=0; i<m; i++ ) { 1303 ourlens[i] -= offlens[i]; 1304 } 1305 ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 1306 A = *newmat; 1307 for ( i=0; i<m; i++ ) { 1308 ourlens[i] += offlens[i]; 1309 } 1310 1311 if (!rank) { 1312 vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 1313 1314 /* read in my part of the matrix numerical values */ 1315 nz = procsnz[0]; 1316 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1317 1318 /* insert into matrix */ 1319 jj = rstart; 1320 smycols = mycols; 1321 svals = vals; 1322 for ( i=0; i<m; i++ ) { 1323 ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1324 smycols += ourlens[i]; 1325 svals += ourlens[i]; 1326 jj++; 1327 } 1328 1329 /* read in other processors and ship out */ 1330 for ( i=1; i<size; i++ ) { 1331 nz = procsnz[i]; 1332 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1333 ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 1334 } 1335 PetscFree(procsnz); 1336 } else { 1337 /* receive numeric values */ 1338 vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 1339 1340 /* receive message of values*/ 1341 ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1342 ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1343 if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 1344 1345 /* insert into matrix */ 1346 jj = rstart; 1347 smycols = mycols; 1348 svals = vals; 1349 for ( i=0; i<m; i++ ) { 1350 ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1351 smycols += ourlens[i]; 1352 svals += ourlens[i]; 1353 jj++; 1354 } 1355 } 1356 PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 1357 1358 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1359 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1360 PetscFunctionReturn(0); 1361 } 1362 1363 1364 1365 1366 1367