18965ea79SLois Curfman McInnes #ifndef lint 2*096963f5SLois Curfman McInnes static char vcid[] = "$Id: mpidense.c,v 1.3 1995/10/22 04:19:52 bsmith Exp curfman $"; 38965ea79SLois Curfman McInnes #endif 48965ea79SLois Curfman McInnes 58965ea79SLois Curfman McInnes #include "mpidense.h" 68965ea79SLois Curfman McInnes #include "vec/vecimpl.h" 78965ea79SLois Curfman McInnes #include "inline/spops.h" 88965ea79SLois Curfman McInnes 939ddd567SLois Curfman McInnes static int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n, 108965ea79SLois Curfman McInnes int *idxn,Scalar *v,InsertMode addv) 118965ea79SLois Curfman McInnes { 1239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1339ddd567SLois Curfman McInnes int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 148965ea79SLois Curfman McInnes 1539ddd567SLois Curfman McInnes if (mdn->insertmode != NOT_SET_VALUES && mdn->insertmode != addv) { 1639ddd567SLois Curfman McInnes SETERRQ(1,"MatSetValues_MPIDense:Cannot mix inserts and adds"); 178965ea79SLois Curfman McInnes } 1839ddd567SLois Curfman McInnes mdn->insertmode = addv; 198965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 2039ddd567SLois Curfman McInnes if (idxm[i] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative row"); 2139ddd567SLois Curfman McInnes if (idxm[i] >= mdn->M) SETERRQ(1,"MatSetValues_MPIDense:Row too large"); 228965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 238965ea79SLois Curfman McInnes row = idxm[i] - rstart; 248965ea79SLois Curfman McInnes for ( j=0; j<n; j++ ) { 2539ddd567SLois Curfman McInnes if (idxn[j] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative column"); 2639ddd567SLois Curfman McInnes if (idxn[j] >= mdn->N) 2739ddd567SLois Curfman McInnes SETERRQ(1,"MatSetValues_MPIDense:Column too large"); 2839ddd567SLois Curfman McInnes ierr = MatSetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j,addv); 2939ddd567SLois Curfman McInnes CHKERRQ(ierr); 308965ea79SLois Curfman McInnes } 318965ea79SLois Curfman McInnes } 328965ea79SLois Curfman McInnes else { 3339ddd567SLois Curfman McInnes ierr = StashValues_Private(&mdn->stash,idxm[i],n,idxn,v+i*n,addv); 3439ddd567SLois Curfman McInnes CHKERRQ(ierr); 358965ea79SLois Curfman McInnes } 368965ea79SLois Curfman McInnes } 378965ea79SLois Curfman McInnes return 0; 388965ea79SLois Curfman McInnes } 398965ea79SLois Curfman McInnes 4039ddd567SLois Curfman McInnes static int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 418965ea79SLois Curfman McInnes { 4239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 438965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 4439ddd567SLois Curfman McInnes int size = mdn->size, *owners = mdn->rowners, rank = mdn->rank; 458965ea79SLois Curfman McInnes int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 4639ddd567SLois Curfman McInnes int tag = mat->tag, *owner,*starts,count,ierr; 478965ea79SLois Curfman McInnes InsertMode addv; 4839ddd567SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 498965ea79SLois Curfman McInnes Scalar *rvalues,*svalues; 508965ea79SLois Curfman McInnes 518965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 5239ddd567SLois Curfman McInnes MPI_Allreduce((void *) &mdn->insertmode,(void *) &addv,1,MPI_INT, 5339ddd567SLois Curfman McInnes MPI_BOR,comm); 5439ddd567SLois Curfman McInnes if (addv == (ADD_VALUES|INSERT_VALUES)) { SETERRQ(1, 5539ddd567SLois Curfman McInnes "MatAssemblyBegin_MPIDense:Cannot mix adds/inserts on different procs"); 568965ea79SLois Curfman McInnes } 5739ddd567SLois Curfman McInnes mdn->insertmode = addv; /* in case this processor had no cache */ 588965ea79SLois Curfman McInnes 598965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 608965ea79SLois Curfman McInnes nprocs = (int *) PETSCMALLOC( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 618965ea79SLois Curfman McInnes PetscZero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 6239ddd567SLois Curfman McInnes owner = (int *) PETSCMALLOC( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 6339ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 6439ddd567SLois Curfman McInnes idx = mdn->stash.idx[i]; 658965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 668965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 678965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; break; 688965ea79SLois Curfman McInnes } 698965ea79SLois Curfman McInnes } 708965ea79SLois Curfman McInnes } 718965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 728965ea79SLois Curfman McInnes 738965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 748965ea79SLois Curfman McInnes work = (int *) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(work); 7539ddd567SLois Curfman McInnes MPI_Allreduce((void *) procs,(void *) work,size,MPI_INT,MPI_SUM,comm); 768965ea79SLois Curfman McInnes nreceives = work[rank]; 7739ddd567SLois Curfman McInnes MPI_Allreduce((void *) nprocs,(void *) work,size,MPI_INT,MPI_MAX,comm); 788965ea79SLois Curfman McInnes nmax = work[rank]; 798965ea79SLois Curfman McInnes PETSCFREE(work); 808965ea79SLois Curfman McInnes 818965ea79SLois Curfman McInnes /* post receives: 828965ea79SLois Curfman McInnes 1) each message will consist of ordered pairs 838965ea79SLois Curfman McInnes (global index,value) we store the global index as a double 848965ea79SLois Curfman McInnes to simplify the message passing. 858965ea79SLois Curfman McInnes 2) since we don't know how long each individual message is we 868965ea79SLois Curfman McInnes allocate the largest needed buffer for each receive. Potentially 878965ea79SLois Curfman McInnes this is a lot of wasted space. 888965ea79SLois Curfman McInnes 898965ea79SLois Curfman McInnes This could be done better. 908965ea79SLois Curfman McInnes */ 918965ea79SLois Curfman McInnes rvalues = (Scalar *) PETSCMALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 928965ea79SLois Curfman McInnes CHKPTRQ(rvalues); 938965ea79SLois Curfman McInnes recv_waits = (MPI_Request *) PETSCMALLOC((nreceives+1)*sizeof(MPI_Request)); 948965ea79SLois Curfman McInnes CHKPTRQ(recv_waits); 958965ea79SLois Curfman McInnes for ( i=0; i<nreceives; i++ ) { 9639ddd567SLois Curfman McInnes MPI_Irecv((void *)(rvalues+3*nmax*i),3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 978965ea79SLois Curfman McInnes comm,recv_waits+i); 988965ea79SLois Curfman McInnes } 998965ea79SLois Curfman McInnes 1008965ea79SLois Curfman McInnes /* do sends: 1018965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 1028965ea79SLois Curfman McInnes the ith processor 1038965ea79SLois Curfman McInnes */ 10439ddd567SLois Curfman McInnes svalues = (Scalar *) PETSCMALLOC( 3*(mdn->stash.n+1)*sizeof(Scalar) ); 10539ddd567SLois Curfman McInnes CHKPTRQ(svalues); 1068965ea79SLois Curfman McInnes send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request)); 1078965ea79SLois Curfman McInnes CHKPTRQ(send_waits); 1088965ea79SLois Curfman McInnes starts = (int *) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(starts); 1098965ea79SLois Curfman McInnes starts[0] = 0; 1108965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 11139ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 11239ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]] = (Scalar) mdn->stash.idx[i]; 11339ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]+1] = (Scalar) mdn->stash.idy[i]; 11439ddd567SLois Curfman McInnes svalues[3*(starts[owner[i]]++)+2] = mdn->stash.array[i]; 1158965ea79SLois Curfman McInnes } 1168965ea79SLois Curfman McInnes PETSCFREE(owner); 1178965ea79SLois Curfman McInnes starts[0] = 0; 1188965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1198965ea79SLois Curfman McInnes count = 0; 1208965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 1218965ea79SLois Curfman McInnes if (procs[i]) { 12239ddd567SLois Curfman McInnes MPI_Isend((void*)(svalues+3*starts[i]),3*nprocs[i],MPIU_SCALAR,i,tag, 1238965ea79SLois Curfman McInnes comm,send_waits+count++); 1248965ea79SLois Curfman McInnes } 1258965ea79SLois Curfman McInnes } 1268965ea79SLois Curfman McInnes PETSCFREE(starts); PETSCFREE(nprocs); 1278965ea79SLois Curfman McInnes 1288965ea79SLois Curfman McInnes /* Free cache space */ 12939ddd567SLois Curfman McInnes ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr); 1308965ea79SLois Curfman McInnes 13139ddd567SLois Curfman McInnes mdn->svalues = svalues; mdn->rvalues = rvalues; 13239ddd567SLois Curfman McInnes mdn->nsends = nsends; mdn->nrecvs = nreceives; 13339ddd567SLois Curfman McInnes mdn->send_waits = send_waits; mdn->recv_waits = recv_waits; 13439ddd567SLois Curfman McInnes mdn->rmax = nmax; 1358965ea79SLois Curfman McInnes 1368965ea79SLois Curfman McInnes return 0; 1378965ea79SLois Curfman McInnes } 13839ddd567SLois Curfman McInnes extern int MatSetUpMultiply_MPIDense(Mat); 1398965ea79SLois Curfman McInnes 14039ddd567SLois Curfman McInnes static int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 1418965ea79SLois Curfman McInnes { 14239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1438965ea79SLois Curfman McInnes MPI_Status *send_status,recv_status; 14439ddd567SLois Curfman McInnes int imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col; 1458965ea79SLois Curfman McInnes Scalar *values,val; 14639ddd567SLois Curfman McInnes InsertMode addv = mdn->insertmode; 1478965ea79SLois Curfman McInnes 1488965ea79SLois Curfman McInnes /* wait on receives */ 1498965ea79SLois Curfman McInnes while (count) { 15039ddd567SLois Curfman McInnes MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status); 1518965ea79SLois Curfman McInnes /* unpack receives into our local space */ 15239ddd567SLois Curfman McInnes values = mdn->rvalues + 3*imdex*mdn->rmax; 1538965ea79SLois Curfman McInnes MPI_Get_count(&recv_status,MPIU_SCALAR,&n); 1548965ea79SLois Curfman McInnes n = n/3; 1558965ea79SLois Curfman McInnes for ( i=0; i<n; i++ ) { 15639ddd567SLois Curfman McInnes row = (int) PETSCREAL(values[3*i]) - mdn->rstart; 1578965ea79SLois Curfman McInnes col = (int) PETSCREAL(values[3*i+1]); 1588965ea79SLois Curfman McInnes val = values[3*i+2]; 15939ddd567SLois Curfman McInnes if (col >= 0 && col < mdn->N) { 16039ddd567SLois Curfman McInnes MatSetValues(mdn->A,1,&row,1,&col,&val,addv); 1618965ea79SLois Curfman McInnes } 16239ddd567SLois Curfman McInnes else {SETERRQ(1,"MatAssemblyEnd_MPIDense:Invalid column");} 1638965ea79SLois Curfman McInnes } 1648965ea79SLois Curfman McInnes count--; 1658965ea79SLois Curfman McInnes } 16639ddd567SLois Curfman McInnes PETSCFREE(mdn->recv_waits); PETSCFREE(mdn->rvalues); 1678965ea79SLois Curfman McInnes 1688965ea79SLois Curfman McInnes /* wait on sends */ 16939ddd567SLois Curfman McInnes if (mdn->nsends) { 17039ddd567SLois Curfman McInnes send_status = (MPI_Status *) PETSCMALLOC( mdn->nsends*sizeof(MPI_Status) ); 1718965ea79SLois Curfman McInnes CHKPTRQ(send_status); 17239ddd567SLois Curfman McInnes MPI_Waitall(mdn->nsends,mdn->send_waits,send_status); 1738965ea79SLois Curfman McInnes PETSCFREE(send_status); 1748965ea79SLois Curfman McInnes } 17539ddd567SLois Curfman McInnes PETSCFREE(mdn->send_waits); PETSCFREE(mdn->svalues); 1768965ea79SLois Curfman McInnes 17739ddd567SLois Curfman McInnes mdn->insertmode = NOT_SET_VALUES; 17839ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr); 17939ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr); 1808965ea79SLois Curfman McInnes 18139ddd567SLois Curfman McInnes if (!mdn->assembled && mode == FINAL_ASSEMBLY) { 18239ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr); 1838965ea79SLois Curfman McInnes } 18439ddd567SLois Curfman McInnes mdn->assembled = 1; 1858965ea79SLois Curfman McInnes return 0; 1868965ea79SLois Curfman McInnes } 1878965ea79SLois Curfman McInnes 18839ddd567SLois Curfman McInnes static int MatZeroEntries_MPIDense(Mat A) 1898965ea79SLois Curfman McInnes { 19039ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 19139ddd567SLois Curfman McInnes return MatZeroEntries(l->A); 1928965ea79SLois Curfman McInnes } 1938965ea79SLois Curfman McInnes 1948965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 1958965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 1968965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 1978965ea79SLois Curfman McInnes aij->A and aij->B directly and not through the MatZeroRows() 1988965ea79SLois Curfman McInnes routine. 1998965ea79SLois Curfman McInnes */ 20039ddd567SLois Curfman McInnes static int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 2018965ea79SLois Curfman McInnes { 20239ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 2038965ea79SLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 2048965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 2058965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2068965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2078965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2088965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2098965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2108965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2118965ea79SLois Curfman McInnes IS istmp; 2128965ea79SLois Curfman McInnes 21339ddd567SLois Curfman McInnes if (!l->assembled) SETERRQ(1,"MatZeroRows_MPIDense:Must assemble matrix"); 2148965ea79SLois Curfman McInnes ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 2158965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 2168965ea79SLois Curfman McInnes 2178965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2188965ea79SLois Curfman McInnes nprocs = (int *) PETSCMALLOC( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 2198965ea79SLois Curfman McInnes PetscZero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 2208965ea79SLois Curfman McInnes owner = (int *) PETSCMALLOC((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 2218965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2228965ea79SLois Curfman McInnes idx = rows[i]; 2238965ea79SLois Curfman McInnes found = 0; 2248965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 2258965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2268965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2278965ea79SLois Curfman McInnes } 2288965ea79SLois Curfman McInnes } 22939ddd567SLois Curfman McInnes if (!found) SETERRQ(1,"MatZeroRows_MPIDense:Index out of range"); 2308965ea79SLois Curfman McInnes } 2318965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 2328965ea79SLois Curfman McInnes 2338965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 2348965ea79SLois Curfman McInnes work = (int *) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(work); 2358965ea79SLois Curfman McInnes MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm); 2368965ea79SLois Curfman McInnes nrecvs = work[rank]; 2378965ea79SLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 2388965ea79SLois Curfman McInnes nmax = work[rank]; 2398965ea79SLois Curfman McInnes PETSCFREE(work); 2408965ea79SLois Curfman McInnes 2418965ea79SLois Curfman McInnes /* post receives: */ 2428965ea79SLois Curfman McInnes rvalues = (int *) PETSCMALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 2438965ea79SLois Curfman McInnes CHKPTRQ(rvalues); 2448965ea79SLois Curfman McInnes recv_waits = (MPI_Request *) PETSCMALLOC((nrecvs+1)*sizeof(MPI_Request)); 2458965ea79SLois Curfman McInnes CHKPTRQ(recv_waits); 2468965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 2478965ea79SLois Curfman McInnes MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i); 2488965ea79SLois Curfman McInnes } 2498965ea79SLois Curfman McInnes 2508965ea79SLois Curfman McInnes /* do sends: 2518965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 2528965ea79SLois Curfman McInnes the ith processor 2538965ea79SLois Curfman McInnes */ 2548965ea79SLois Curfman McInnes svalues = (int *) PETSCMALLOC( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 2558965ea79SLois Curfman McInnes send_waits = (MPI_Request *) PETSCMALLOC( (nsends+1)*sizeof(MPI_Request)); 2568965ea79SLois Curfman McInnes CHKPTRQ(send_waits); 2578965ea79SLois Curfman McInnes starts = (int *) PETSCMALLOC( (size+1)*sizeof(int) ); CHKPTRQ(starts); 2588965ea79SLois Curfman McInnes starts[0] = 0; 2598965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2608965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2618965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 2628965ea79SLois Curfman McInnes } 2638965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 2648965ea79SLois Curfman McInnes 2658965ea79SLois Curfman McInnes starts[0] = 0; 2668965ea79SLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2678965ea79SLois Curfman McInnes count = 0; 2688965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 2698965ea79SLois Curfman McInnes if (procs[i]) { 2708965ea79SLois Curfman McInnes MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++); 2718965ea79SLois Curfman McInnes } 2728965ea79SLois Curfman McInnes } 2738965ea79SLois Curfman McInnes PETSCFREE(starts); 2748965ea79SLois Curfman McInnes 2758965ea79SLois Curfman McInnes base = owners[rank]; 2768965ea79SLois Curfman McInnes 2778965ea79SLois Curfman McInnes /* wait on receives */ 2788965ea79SLois Curfman McInnes lens = (int *) PETSCMALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 2798965ea79SLois Curfman McInnes source = lens + nrecvs; 2808965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 2818965ea79SLois Curfman McInnes while (count) { 2828965ea79SLois Curfman McInnes MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 2838965ea79SLois Curfman McInnes /* unpack receives into our local space */ 2848965ea79SLois Curfman McInnes MPI_Get_count(&recv_status,MPI_INT,&n); 2858965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 2868965ea79SLois Curfman McInnes lens[imdex] = n; 2878965ea79SLois Curfman McInnes slen += n; 2888965ea79SLois Curfman McInnes count--; 2898965ea79SLois Curfman McInnes } 2908965ea79SLois Curfman McInnes PETSCFREE(recv_waits); 2918965ea79SLois Curfman McInnes 2928965ea79SLois Curfman McInnes /* move the data into the send scatter */ 2938965ea79SLois Curfman McInnes lrows = (int *) PETSCMALLOC( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 2948965ea79SLois Curfman McInnes count = 0; 2958965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 2968965ea79SLois Curfman McInnes values = rvalues + i*nmax; 2978965ea79SLois Curfman McInnes for ( j=0; j<lens[i]; j++ ) { 2988965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 2998965ea79SLois Curfman McInnes } 3008965ea79SLois Curfman McInnes } 3018965ea79SLois Curfman McInnes PETSCFREE(rvalues); PETSCFREE(lens); 3028965ea79SLois Curfman McInnes PETSCFREE(owner); PETSCFREE(nprocs); 3038965ea79SLois Curfman McInnes 3048965ea79SLois Curfman McInnes /* actually zap the local rows */ 3058965ea79SLois Curfman McInnes ierr = ISCreateSeq(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 3068965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 3078965ea79SLois Curfman McInnes PETSCFREE(lrows); 3088965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 3098965ea79SLois Curfman McInnes ierr = ISDestroy(istmp); CHKERRQ(ierr); 3108965ea79SLois Curfman McInnes 3118965ea79SLois Curfman McInnes /* wait on sends */ 3128965ea79SLois Curfman McInnes if (nsends) { 3138965ea79SLois Curfman McInnes send_status = (MPI_Status *) PETSCMALLOC(nsends*sizeof(MPI_Status)); 3148965ea79SLois Curfman McInnes CHKPTRQ(send_status); 3158965ea79SLois Curfman McInnes MPI_Waitall(nsends,send_waits,send_status); 3168965ea79SLois Curfman McInnes PETSCFREE(send_status); 3178965ea79SLois Curfman McInnes } 3188965ea79SLois Curfman McInnes PETSCFREE(send_waits); PETSCFREE(svalues); 3198965ea79SLois Curfman McInnes 3208965ea79SLois Curfman McInnes return 0; 3218965ea79SLois Curfman McInnes } 3228965ea79SLois Curfman McInnes 32339ddd567SLois Curfman McInnes static int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3248965ea79SLois Curfman McInnes { 32539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3268965ea79SLois Curfman McInnes int ierr; 32739ddd567SLois Curfman McInnes if (!mdn->assembled) 32839ddd567SLois Curfman McInnes SETERRQ(1,"MatMult_MPIDense:Must assemble matrix first"); 32939ddd567SLois Curfman McInnes ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 33039ddd567SLois Curfman McInnes CHKERRQ(ierr); 33139ddd567SLois Curfman McInnes ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 33239ddd567SLois Curfman McInnes CHKERRQ(ierr); 33339ddd567SLois Curfman McInnes ierr = MatMult(mdn->A,mdn->lvec,yy); CHKERRQ(ierr); 3348965ea79SLois Curfman McInnes return 0; 3358965ea79SLois Curfman McInnes } 3368965ea79SLois Curfman McInnes 33739ddd567SLois Curfman McInnes static int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3388965ea79SLois Curfman McInnes { 33939ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3408965ea79SLois Curfman McInnes int ierr; 34139ddd567SLois Curfman McInnes if (!mdn->assembled) 34239ddd567SLois Curfman McInnes SETERRQ(1,"MatMultAdd_MPIDense:Must assemble matrix first"); 34339ddd567SLois Curfman McInnes ierr = VecScatterBegin(xx,mdn->lvec,ADD_VALUES,SCATTER_ALL,mdn->Mvctx); 34439ddd567SLois Curfman McInnes CHKERRQ(ierr); 34539ddd567SLois Curfman McInnes ierr = VecScatterEnd(xx,mdn->lvec,ADD_VALUES,SCATTER_ALL,mdn->Mvctx); 34639ddd567SLois Curfman McInnes CHKERRQ(ierr); 34739ddd567SLois Curfman McInnes ierr = MatMultAdd(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr); 3488965ea79SLois Curfman McInnes return 0; 3498965ea79SLois Curfman McInnes } 3508965ea79SLois Curfman McInnes 351*096963f5SLois Curfman McInnes static int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 352*096963f5SLois Curfman McInnes { 353*096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 354*096963f5SLois Curfman McInnes int ierr; 355*096963f5SLois Curfman McInnes 356*096963f5SLois Curfman McInnes if (!a->assembled) SETERRQ(1,"MatMulTrans_MPIDense:must assemble matrix"); 357*096963f5SLois Curfman McInnes ierr = MatMultTrans(a->A,xx,a->lvec); CHKERRQ(ierr); 358*096963f5SLois Curfman McInnes ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES, 359*096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 360*096963f5SLois Curfman McInnes ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES, 361*096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 362*096963f5SLois Curfman McInnes return 0; 363*096963f5SLois Curfman McInnes } 364*096963f5SLois Curfman McInnes 365*096963f5SLois Curfman McInnes static int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 366*096963f5SLois Curfman McInnes { 367*096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 368*096963f5SLois Curfman McInnes int ierr; 369*096963f5SLois Curfman McInnes 370*096963f5SLois Curfman McInnes if (!a->assembled) SETERRQ(1,"MatMulTransAdd_MPIDense:must assemble matrix"); 371*096963f5SLois Curfman McInnes ierr = MatMultTransAdd(a->A,xx,yy,a->lvec); CHKERRQ(ierr); 372*096963f5SLois Curfman McInnes /* send it on its way */ 373*096963f5SLois Curfman McInnes ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES, 374*096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 375*096963f5SLois Curfman McInnes ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES, 376*096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 377*096963f5SLois Curfman McInnes return 0; 378*096963f5SLois Curfman McInnes } 379*096963f5SLois Curfman McInnes 38039ddd567SLois Curfman McInnes static int MatGetDiagonal_MPIDense(Mat A,Vec v) 3818965ea79SLois Curfman McInnes { 38239ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 383*096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 384*096963f5SLois Curfman McInnes int ierr, i, n, m = a->m, radd; 385*096963f5SLois Curfman McInnes Scalar *x; 38639ddd567SLois Curfman McInnes if (!a->assembled) SETERRQ(1,"MatGetDiag_MPIDense:must assemble matrix"); 387*096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x); CHKERRQ(ierr); 388*096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n); CHKERRQ(ierr); 389*096963f5SLois Curfman McInnes if (n != a->M) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 390*096963f5SLois Curfman McInnes radd = a->rstart*m*m; 391*096963f5SLois Curfman McInnes for ( i=0; i<m; i++ ) { 392*096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 393*096963f5SLois Curfman McInnes } 394*096963f5SLois Curfman McInnes return 0; 3958965ea79SLois Curfman McInnes } 3968965ea79SLois Curfman McInnes 39739ddd567SLois Curfman McInnes static int MatDestroy_MPIDense(PetscObject obj) 3988965ea79SLois Curfman McInnes { 3998965ea79SLois Curfman McInnes Mat mat = (Mat) obj; 40039ddd567SLois Curfman McInnes Mat_MPIDense *aij = (Mat_MPIDense *) mat->data; 4018965ea79SLois Curfman McInnes int ierr; 4028965ea79SLois Curfman McInnes #if defined(PETSC_LOG) 4038965ea79SLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N); 4048965ea79SLois Curfman McInnes #endif 4058965ea79SLois Curfman McInnes PETSCFREE(aij->rowners); 4068965ea79SLois Curfman McInnes ierr = MatDestroy(aij->A); CHKERRQ(ierr); 4078965ea79SLois Curfman McInnes if (aij->lvec) VecDestroy(aij->lvec); 40855659b69SBarry Smith if (aij->Mvctx) VecScatterDestroy(aij->Mvctx); 4098965ea79SLois Curfman McInnes PETSCFREE(aij); 4108965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4118965ea79SLois Curfman McInnes PETSCHEADERDESTROY(mat); 4128965ea79SLois Curfman McInnes return 0; 4138965ea79SLois Curfman McInnes } 41439ddd567SLois Curfman McInnes 4158965ea79SLois Curfman McInnes #include "pinclude/pviewer.h" 4168965ea79SLois Curfman McInnes 41739ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 4188965ea79SLois Curfman McInnes { 41939ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4208965ea79SLois Curfman McInnes int ierr; 42139ddd567SLois Curfman McInnes if (mdn->size == 1) { 42239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 4238965ea79SLois Curfman McInnes } 42439ddd567SLois Curfman McInnes else SETERRQ(1,"MatView_MPIDense_Binary:Only uniprocessor output supported"); 4258965ea79SLois Curfman McInnes return 0; 4268965ea79SLois Curfman McInnes } 4278965ea79SLois Curfman McInnes 42839ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 4298965ea79SLois Curfman McInnes { 43039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 43139ddd567SLois Curfman McInnes int ierr, format; 4328965ea79SLois Curfman McInnes PetscObject vobj = (PetscObject) viewer; 4338965ea79SLois Curfman McInnes FILE *fd; 4348965ea79SLois Curfman McInnes 4358965ea79SLois Curfman McInnes if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 4368965ea79SLois Curfman McInnes ierr = ViewerFileGetFormat_Private(viewer,&format); 4378965ea79SLois Curfman McInnes if (format == FILE_FORMAT_INFO) { 438*096963f5SLois Curfman McInnes int nz, nzalloc, mem, rank; 439*096963f5SLois Curfman McInnes MPI_Comm_rank(mat->comm,&rank); 440*096963f5SLois Curfman McInnes ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 441*096963f5SLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&nz,&nzalloc,&mem); 442*096963f5SLois Curfman McInnes MPIU_Seq_begin(mat->comm,1); 443*096963f5SLois Curfman McInnes fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d \n", 444*096963f5SLois Curfman McInnes rank,mdn->m,nz,nzalloc,mem); 445*096963f5SLois Curfman McInnes fflush(fd); 446*096963f5SLois Curfman McInnes MPIU_Seq_end(mat->comm,1); 447*096963f5SLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer); 448*096963f5SLois Curfman McInnes return 0; 4498965ea79SLois Curfman McInnes } 4508965ea79SLois Curfman McInnes } 4518965ea79SLois Curfman McInnes 4528965ea79SLois Curfman McInnes if (vobj->type == ASCII_FILE_VIEWER) { 4538965ea79SLois Curfman McInnes ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 4548965ea79SLois Curfman McInnes MPIU_Seq_begin(mat->comm,1); 45539ddd567SLois Curfman McInnes fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n", 45639ddd567SLois Curfman McInnes mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n); 45739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 4588965ea79SLois Curfman McInnes fflush(fd); 4598965ea79SLois Curfman McInnes MPIU_Seq_end(mat->comm,1); 4608965ea79SLois Curfman McInnes } 4618965ea79SLois Curfman McInnes else { 46239ddd567SLois Curfman McInnes int size = mdn->size, rank = mdn->rank; 4638965ea79SLois Curfman McInnes if (size == 1) { 46439ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 4658965ea79SLois Curfman McInnes } 4668965ea79SLois Curfman McInnes else { 4678965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 4688965ea79SLois Curfman McInnes Mat A; 46939ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 47039ddd567SLois Curfman McInnes Scalar *vals; 47139ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 4728965ea79SLois Curfman McInnes 4738965ea79SLois Curfman McInnes if (!rank) { 47439ddd567SLois Curfman McInnes ierr = MatCreateMPIDense(mat->comm,M,M,N,N,&A); CHKERRQ(ierr); 4758965ea79SLois Curfman McInnes } 4768965ea79SLois Curfman McInnes else { 47739ddd567SLois Curfman McInnes ierr = MatCreateMPIDense(mat->comm,0,M,N,N,&A); CHKERRQ(ierr); 4788965ea79SLois Curfman McInnes } 4798965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 4808965ea79SLois Curfman McInnes 48139ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 48239ddd567SLois Curfman McInnes but it's quick for now */ 48339ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 4848965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 48539ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 48639ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 48739ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 48839ddd567SLois Curfman McInnes row++; 4898965ea79SLois Curfman McInnes } 4908965ea79SLois Curfman McInnes 4918965ea79SLois Curfman McInnes ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 4928965ea79SLois Curfman McInnes ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 4938965ea79SLois Curfman McInnes if (!rank) { 49439ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr); 4958965ea79SLois Curfman McInnes } 4968965ea79SLois Curfman McInnes ierr = MatDestroy(A); CHKERRQ(ierr); 4978965ea79SLois Curfman McInnes } 4988965ea79SLois Curfman McInnes } 4998965ea79SLois Curfman McInnes return 0; 5008965ea79SLois Curfman McInnes } 5018965ea79SLois Curfman McInnes 50239ddd567SLois Curfman McInnes static int MatView_MPIDense(PetscObject obj,Viewer viewer) 5038965ea79SLois Curfman McInnes { 5048965ea79SLois Curfman McInnes Mat mat = (Mat) obj; 50539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 5068965ea79SLois Curfman McInnes PetscObject vobj = (PetscObject) viewer; 50739ddd567SLois Curfman McInnes int ierr; 5088965ea79SLois Curfman McInnes 50939ddd567SLois Curfman McInnes if (!mdn->assembled) SETERRQ(1,"MatView_MPIDense:must assemble matrix"); 5108965ea79SLois Curfman McInnes if (!viewer) { 5118965ea79SLois Curfman McInnes viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 5128965ea79SLois Curfman McInnes } 51339ddd567SLois Curfman McInnes if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILE_VIEWER) { 51439ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 5158965ea79SLois Curfman McInnes } 5168965ea79SLois Curfman McInnes else if (vobj->cookie == VIEWER_COOKIE && vobj->type == ASCII_FILES_VIEWER) { 51739ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 5188965ea79SLois Curfman McInnes } 5198965ea79SLois Curfman McInnes else if (vobj->type == BINARY_FILE_VIEWER) { 52039ddd567SLois Curfman McInnes return MatView_MPIDense_Binary(mat,viewer); 5218965ea79SLois Curfman McInnes } 5228965ea79SLois Curfman McInnes return 0; 5238965ea79SLois Curfman McInnes } 5248965ea79SLois Curfman McInnes 52539ddd567SLois Curfman McInnes static int MatGetInfo_MPIDense(Mat matin,MatInfoType flag,int *nz, 5268965ea79SLois Curfman McInnes int *nzalloc,int *mem) 5278965ea79SLois Curfman McInnes { 52839ddd567SLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) matin->data; 52939ddd567SLois Curfman McInnes Mat A = mat->A; 53039ddd567SLois Curfman McInnes int ierr, isend[3], irecv[3]; 5318965ea79SLois Curfman McInnes 53239ddd567SLois Curfman McInnes ierr = MatGetInfo(A,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr); 5338965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 5348965ea79SLois Curfman McInnes *nz = isend[0]; *nzalloc = isend[1]; *mem = isend[2]; 5358965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 5368965ea79SLois Curfman McInnes MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_MAX,matin->comm); 5378965ea79SLois Curfman McInnes *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2]; 5388965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 5398965ea79SLois Curfman McInnes MPI_Allreduce( isend, irecv,3,MPI_INT,MPI_SUM,matin->comm); 5408965ea79SLois Curfman McInnes *nz = irecv[0]; *nzalloc = irecv[1]; *mem = irecv[2]; 5418965ea79SLois Curfman McInnes } 5428965ea79SLois Curfman McInnes return 0; 5438965ea79SLois Curfman McInnes } 5448965ea79SLois Curfman McInnes 54539ddd567SLois Curfman McInnes static int MatSetOption_MPIDense(Mat A,MatOption op) 5468965ea79SLois Curfman McInnes { 54739ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 5488965ea79SLois Curfman McInnes 5498965ea79SLois Curfman McInnes if (op == NO_NEW_NONZERO_LOCATIONS || 5508965ea79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 5518965ea79SLois Curfman McInnes op == COLUMNS_SORTED || 5528965ea79SLois Curfman McInnes op == ROW_ORIENTED) { 5538965ea79SLois Curfman McInnes MatSetOption(a->A,op); 5548965ea79SLois Curfman McInnes } 5558965ea79SLois Curfman McInnes else if (op == ROWS_SORTED || 5568965ea79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 5578965ea79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 5588965ea79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 55939ddd567SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_MPIDense:Option ignored\n"); 5608965ea79SLois Curfman McInnes else if (op == COLUMN_ORIENTED) 56139ddd567SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:COLUMN_ORIENTED");} 5628965ea79SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 56339ddd567SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:NO_NEW_DIAGONALS");} 5648965ea79SLois Curfman McInnes else 56539ddd567SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:unknown option");} 5668965ea79SLois Curfman McInnes return 0; 5678965ea79SLois Curfman McInnes } 5688965ea79SLois Curfman McInnes 56939ddd567SLois Curfman McInnes static int MatGetSize_MPIDense(Mat matin,int *m,int *n) 5708965ea79SLois Curfman McInnes { 57139ddd567SLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) matin->data; 5728965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 5738965ea79SLois Curfman McInnes return 0; 5748965ea79SLois Curfman McInnes } 5758965ea79SLois Curfman McInnes 57639ddd567SLois Curfman McInnes static int MatGetLocalSize_MPIDense(Mat matin,int *m,int *n) 5778965ea79SLois Curfman McInnes { 57839ddd567SLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) matin->data; 5798965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 5808965ea79SLois Curfman McInnes return 0; 5818965ea79SLois Curfman McInnes } 5828965ea79SLois Curfman McInnes 58339ddd567SLois Curfman McInnes static int MatGetOwnershipRange_MPIDense(Mat matin,int *m,int *n) 5848965ea79SLois Curfman McInnes { 58539ddd567SLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) matin->data; 5868965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 5878965ea79SLois Curfman McInnes return 0; 5888965ea79SLois Curfman McInnes } 5898965ea79SLois Curfman McInnes 59039ddd567SLois Curfman McInnes static int MatGetRow_MPIDense(Mat matin,int row,int *nz,int **idx,Scalar **v) 5918965ea79SLois Curfman McInnes { 59239ddd567SLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) matin->data; 59339ddd567SLois Curfman McInnes int lrow, rstart = mat->rstart, rend = mat->rend; 5948965ea79SLois Curfman McInnes 59539ddd567SLois Curfman McInnes if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIDense:only local rows") 5968965ea79SLois Curfman McInnes lrow = row - rstart; 59739ddd567SLois Curfman McInnes return MatGetRow(mat->A,lrow,nz,idx,v); 5988965ea79SLois Curfman McInnes } 5998965ea79SLois Curfman McInnes 60039ddd567SLois Curfman McInnes static int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 6018965ea79SLois Curfman McInnes { 6028965ea79SLois Curfman McInnes if (idx) PETSCFREE(*idx); 6038965ea79SLois Curfman McInnes if (v) PETSCFREE(*v); 6048965ea79SLois Curfman McInnes return 0; 6058965ea79SLois Curfman McInnes } 6068965ea79SLois Curfman McInnes 607*096963f5SLois Curfman McInnes static int MatNorm_MPIDense(Mat mat,MatNormType type,double *norm) 608*096963f5SLois Curfman McInnes { 609*096963f5SLois Curfman McInnes *norm = 0.0; 610*096963f5SLois Curfman McInnes SETERRQ(1,"MatNorm_MPIDense:Not finished"); 611*096963f5SLois Curfman McInnes return 0; 612*096963f5SLois Curfman McInnes } 613*096963f5SLois Curfman McInnes 61455659b69SBarry Smith static int MatCopyPrivate_MPIDense(Mat,Mat *,int); 6158965ea79SLois Curfman McInnes 6168965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 61739ddd567SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIDense, 61839ddd567SLois Curfman McInnes MatGetRow_MPIDense,MatRestoreRow_MPIDense, 61939ddd567SLois Curfman McInnes MatMult_MPIDense,MatMultAdd_MPIDense, 620*096963f5SLois Curfman McInnes MatMultTrans_MPIDense,MatMultTransAdd_MPIDense, 62139ddd567SLois Curfman McInnes 0, 0, 62239ddd567SLois Curfman McInnes 0, 0, 0, 62339ddd567SLois Curfman McInnes 0, 0, 0, 62439ddd567SLois Curfman McInnes MatGetInfo_MPIDense, 0, 625*096963f5SLois Curfman McInnes MatGetDiagonal_MPIDense, 0, MatNorm_MPIDense, 62639ddd567SLois Curfman McInnes MatAssemblyBegin_MPIDense,MatAssemblyEnd_MPIDense, 6278965ea79SLois Curfman McInnes 0, 62839ddd567SLois Curfman McInnes MatSetOption_MPIDense,MatZeroEntries_MPIDense,MatZeroRows_MPIDense, 62939ddd567SLois Curfman McInnes 0, 63039ddd567SLois Curfman McInnes 0, 0, 0, 0, 63139ddd567SLois Curfman McInnes MatGetSize_MPIDense,MatGetLocalSize_MPIDense, 63239ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 63339ddd567SLois Curfman McInnes 0, 0, 63439ddd567SLois Curfman McInnes 0, 0, 0, 0, 0, MatCopyPrivate_MPIDense}; 6358965ea79SLois Curfman McInnes 6368965ea79SLois Curfman McInnes /*@C 63739ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 6388965ea79SLois Curfman McInnes 6398965ea79SLois Curfman McInnes Input Parameters: 6408965ea79SLois Curfman McInnes . comm - MPI communicator 6418965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 6428965ea79SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated 6438965ea79SLois Curfman McInnes if N is given) 6448965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 6458965ea79SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated 6468965ea79SLois Curfman McInnes if n is given) 6478965ea79SLois Curfman McInnes 6488965ea79SLois Curfman McInnes Output Parameter: 6498965ea79SLois Curfman McInnes . newmat - the matrix 6508965ea79SLois Curfman McInnes 6518965ea79SLois Curfman McInnes Notes: 65239ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 65339ddd567SLois Curfman McInnes storage by columns. 6548965ea79SLois Curfman McInnes 6558965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 6568965ea79SLois Curfman McInnes (possibly both). 6578965ea79SLois Curfman McInnes 65839ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 6598965ea79SLois Curfman McInnes 66039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 6618965ea79SLois Curfman McInnes @*/ 66239ddd567SLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Mat *newmat) 6638965ea79SLois Curfman McInnes { 6648965ea79SLois Curfman McInnes Mat mat; 66539ddd567SLois Curfman McInnes Mat_MPIDense *a; 66639ddd567SLois Curfman McInnes int ierr, i; 6678965ea79SLois Curfman McInnes 6688965ea79SLois Curfman McInnes *newmat = 0; 66939ddd567SLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIDENSE,comm); 6708965ea79SLois Curfman McInnes PLogObjectCreate(mat); 67139ddd567SLois Curfman McInnes mat->data = (void *) (a = PETSCNEW(Mat_MPIDense)); CHKPTRQ(a); 6728965ea79SLois Curfman McInnes PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 67339ddd567SLois Curfman McInnes mat->destroy = MatDestroy_MPIDense; 67439ddd567SLois Curfman McInnes mat->view = MatView_MPIDense; 6758965ea79SLois Curfman McInnes mat->factor = 0; 6768965ea79SLois Curfman McInnes 6778965ea79SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 6788965ea79SLois Curfman McInnes MPI_Comm_rank(comm,&a->rank); 6798965ea79SLois Curfman McInnes MPI_Comm_size(comm,&a->size); 6808965ea79SLois Curfman McInnes 68139ddd567SLois Curfman McInnes if (M == PETSC_DECIDE) MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm); 6828965ea79SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);} 68339ddd567SLois Curfman McInnes 68439ddd567SLois Curfman McInnes /* each row stores all columns */ 68539ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 68639ddd567SLois Curfman McInnes if (n == PETSC_DECIDE) n = N; 68739ddd567SLois Curfman McInnes if (n != N) SETERRQ(1,"MatCreateMPIDense:For now, only n=N is supported"); 6888965ea79SLois Curfman McInnes a->N = N; 6898965ea79SLois Curfman McInnes a->M = M; 69039ddd567SLois Curfman McInnes a->m = m; 69139ddd567SLois Curfman McInnes a->n = n; 6928965ea79SLois Curfman McInnes 6938965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 69439ddd567SLois Curfman McInnes a->rowners = (int *) PETSCMALLOC((a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 69539ddd567SLois Curfman McInnes PLogObjectMemory(mat,(a->size+2)*sizeof(int)+sizeof(struct _Mat)+ 69639ddd567SLois Curfman McInnes sizeof(Mat_MPIDense)); 6978965ea79SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm); 6988965ea79SLois Curfman McInnes a->rowners[0] = 0; 6998965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 7008965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 7018965ea79SLois Curfman McInnes } 7028965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 7038965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 7048965ea79SLois Curfman McInnes 70539ddd567SLois Curfman McInnes ierr = MatCreateSeqDense(MPI_COMM_SELF,m,N,&a->A); CHKERRQ(ierr); 7068965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 7078965ea79SLois Curfman McInnes 7088965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 7098965ea79SLois Curfman McInnes ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr); 7108965ea79SLois Curfman McInnes 7118965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 7128965ea79SLois Curfman McInnes a->lvec = 0; 7138965ea79SLois Curfman McInnes a->Mvctx = 0; 7148965ea79SLois Curfman McInnes a->assembled = 0; 7158965ea79SLois Curfman McInnes 7168965ea79SLois Curfman McInnes *newmat = mat; 7178965ea79SLois Curfman McInnes return 0; 7188965ea79SLois Curfman McInnes } 7198965ea79SLois Curfman McInnes 72055659b69SBarry Smith static int MatCopyPrivate_MPIDense(Mat matin,Mat *newmat,int cpvalues) 7218965ea79SLois Curfman McInnes { 7228965ea79SLois Curfman McInnes Mat mat; 72339ddd567SLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) matin->data; 72439ddd567SLois Curfman McInnes int ierr; 7258965ea79SLois Curfman McInnes 72639ddd567SLois Curfman McInnes if (!oldmat->assembled) SETERRQ(1,"MatCopyPrivate_MPIDense:Must assemble matrix"); 7278965ea79SLois Curfman McInnes *newmat = 0; 72839ddd567SLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIDENSE,matin->comm); 7298965ea79SLois Curfman McInnes PLogObjectCreate(mat); 73039ddd567SLois Curfman McInnes mat->data = (void *) (a = PETSCNEW(Mat_MPIDense)); CHKPTRQ(a); 7318965ea79SLois Curfman McInnes PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 73239ddd567SLois Curfman McInnes mat->destroy = MatDestroy_MPIDense; 73339ddd567SLois Curfman McInnes mat->view = MatView_MPIDense; 7348965ea79SLois Curfman McInnes mat->factor = matin->factor; 7358965ea79SLois Curfman McInnes 7368965ea79SLois Curfman McInnes a->m = oldmat->m; 7378965ea79SLois Curfman McInnes a->n = oldmat->n; 7388965ea79SLois Curfman McInnes a->M = oldmat->M; 7398965ea79SLois Curfman McInnes a->N = oldmat->N; 7408965ea79SLois Curfman McInnes 7418965ea79SLois Curfman McInnes a->assembled = 1; 7428965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 7438965ea79SLois Curfman McInnes a->rend = oldmat->rend; 7448965ea79SLois Curfman McInnes a->size = oldmat->size; 7458965ea79SLois Curfman McInnes a->rank = oldmat->rank; 7468965ea79SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 7478965ea79SLois Curfman McInnes 7488965ea79SLois Curfman McInnes a->rowners = (int *) PETSCMALLOC((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners); 74939ddd567SLois Curfman McInnes PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIDense)); 7508965ea79SLois Curfman McInnes PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int)); 7518965ea79SLois Curfman McInnes ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 7528965ea79SLois Curfman McInnes 7538965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 7548965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 75555659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 7568965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 7578965ea79SLois Curfman McInnes ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 7588965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 7598965ea79SLois Curfman McInnes *newmat = mat; 7608965ea79SLois Curfman McInnes return 0; 7618965ea79SLois Curfman McInnes } 7628965ea79SLois Curfman McInnes 7638965ea79SLois Curfman McInnes #include "sysio.h" 7648965ea79SLois Curfman McInnes 76539ddd567SLois Curfman McInnes int MatLoad_MPIDense(Viewer bview,MatType type,Mat *newmat) 7668965ea79SLois Curfman McInnes { 7678965ea79SLois Curfman McInnes Mat A; 7688965ea79SLois Curfman McInnes int i, nz, ierr, j,rstart, rend, fd; 7698965ea79SLois Curfman McInnes Scalar *vals,*svals; 7708965ea79SLois Curfman McInnes PetscObject vobj = (PetscObject) bview; 7718965ea79SLois Curfman McInnes MPI_Comm comm = vobj->comm; 7728965ea79SLois Curfman McInnes MPI_Status status; 7738965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 7748965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 7758965ea79SLois Curfman McInnes int tag = ((PetscObject)bview)->tag; 7768965ea79SLois Curfman McInnes 7778965ea79SLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 7788965ea79SLois Curfman McInnes if (!rank) { 7798965ea79SLois Curfman McInnes ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 7808965ea79SLois Curfman McInnes ierr = SYRead(fd,(char *)header,4,SYINT); CHKERRQ(ierr); 78139ddd567SLois Curfman McInnes if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:not matrix object"); 7828965ea79SLois Curfman McInnes } 7838965ea79SLois Curfman McInnes 7848965ea79SLois Curfman McInnes MPI_Bcast(header+1,3,MPI_INT,0,comm); 7858965ea79SLois Curfman McInnes M = header[1]; N = header[2]; 7868965ea79SLois Curfman McInnes /* determine ownership of all rows */ 7878965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 7888965ea79SLois Curfman McInnes rowners = (int *) PETSCMALLOC((size+2)*sizeof(int)); CHKPTRQ(rowners); 7898965ea79SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 7908965ea79SLois Curfman McInnes rowners[0] = 0; 7918965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 7928965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 7938965ea79SLois Curfman McInnes } 7948965ea79SLois Curfman McInnes rstart = rowners[rank]; 7958965ea79SLois Curfman McInnes rend = rowners[rank+1]; 7968965ea79SLois Curfman McInnes 7978965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 7988965ea79SLois Curfman McInnes ourlens = (int*) PETSCMALLOC( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 7998965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 8008965ea79SLois Curfman McInnes if (!rank) { 8018965ea79SLois Curfman McInnes rowlengths = (int*) PETSCMALLOC( M*sizeof(int) ); CHKPTRQ(rowlengths); 8028965ea79SLois Curfman McInnes ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 8038965ea79SLois Curfman McInnes sndcounts = (int*) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(sndcounts); 8048965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 8058965ea79SLois Curfman McInnes MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm); 8068965ea79SLois Curfman McInnes PETSCFREE(sndcounts); 8078965ea79SLois Curfman McInnes } 8088965ea79SLois Curfman McInnes else { 8098965ea79SLois Curfman McInnes MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm); 8108965ea79SLois Curfman McInnes } 8118965ea79SLois Curfman McInnes 8128965ea79SLois Curfman McInnes if (!rank) { 8138965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 8148965ea79SLois Curfman McInnes procsnz = (int*) PETSCMALLOC( size*sizeof(int) ); CHKPTRQ(procsnz); 8158965ea79SLois Curfman McInnes PetscZero(procsnz,size*sizeof(int)); 8168965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 8178965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 8188965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 8198965ea79SLois Curfman McInnes } 8208965ea79SLois Curfman McInnes } 8218965ea79SLois Curfman McInnes PETSCFREE(rowlengths); 8228965ea79SLois Curfman McInnes 8238965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 8248965ea79SLois Curfman McInnes maxnz = 0; 8258965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 8268965ea79SLois Curfman McInnes maxnz = PETSCMAX(maxnz,procsnz[i]); 8278965ea79SLois Curfman McInnes } 8288965ea79SLois Curfman McInnes cols = (int *) PETSCMALLOC( maxnz*sizeof(int) ); CHKPTRQ(cols); 8298965ea79SLois Curfman McInnes 8308965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 8318965ea79SLois Curfman McInnes nz = procsnz[0]; 8328965ea79SLois Curfman McInnes mycols = (int *) PETSCMALLOC( nz*sizeof(int) ); CHKPTRQ(mycols); 8338965ea79SLois Curfman McInnes ierr = SYRead(fd,mycols,nz,SYINT); CHKERRQ(ierr); 8348965ea79SLois Curfman McInnes 8358965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 8368965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 8378965ea79SLois Curfman McInnes nz = procsnz[i]; 8388965ea79SLois Curfman McInnes ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr); 8398965ea79SLois Curfman McInnes MPI_Send(cols,nz,MPI_INT,i,tag,comm); 8408965ea79SLois Curfman McInnes } 8418965ea79SLois Curfman McInnes PETSCFREE(cols); 8428965ea79SLois Curfman McInnes } 8438965ea79SLois Curfman McInnes else { 8448965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 8458965ea79SLois Curfman McInnes nz = 0; 8468965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 8478965ea79SLois Curfman McInnes nz += ourlens[i]; 8488965ea79SLois Curfman McInnes } 8498965ea79SLois Curfman McInnes mycols = (int*) PETSCMALLOC( nz*sizeof(int) ); CHKPTRQ(mycols); 8508965ea79SLois Curfman McInnes 8518965ea79SLois Curfman McInnes /* receive message of column indices*/ 8528965ea79SLois Curfman McInnes MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status); 8538965ea79SLois Curfman McInnes MPI_Get_count(&status,MPI_INT,&maxnz); 85439ddd567SLois Curfman McInnes if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file"); 8558965ea79SLois Curfman McInnes } 8568965ea79SLois Curfman McInnes 8578965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 8588965ea79SLois Curfman McInnes PetscZero(offlens,m*sizeof(int)); 8598965ea79SLois Curfman McInnes jj = 0; 8608965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 8618965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 8628965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 8638965ea79SLois Curfman McInnes jj++; 8648965ea79SLois Curfman McInnes } 8658965ea79SLois Curfman McInnes } 8668965ea79SLois Curfman McInnes 8678965ea79SLois Curfman McInnes /* create our matrix */ 8688965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 8698965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 8708965ea79SLois Curfman McInnes } 87139ddd567SLois Curfman McInnes if (type == MATMPIDENSE) { 87239ddd567SLois Curfman McInnes ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 8738965ea79SLois Curfman McInnes } 8748965ea79SLois Curfman McInnes A = *newmat; 8758965ea79SLois Curfman McInnes MatSetOption(A,COLUMNS_SORTED); 8768965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 8778965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 8788965ea79SLois Curfman McInnes } 8798965ea79SLois Curfman McInnes 8808965ea79SLois Curfman McInnes if (!rank) { 8818965ea79SLois Curfman McInnes vals = (Scalar *) PETSCMALLOC( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 8828965ea79SLois Curfman McInnes 8838965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 8848965ea79SLois Curfman McInnes nz = procsnz[0]; 8858965ea79SLois Curfman McInnes ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 8868965ea79SLois Curfman McInnes 8878965ea79SLois Curfman McInnes /* insert into matrix */ 8888965ea79SLois Curfman McInnes jj = rstart; 8898965ea79SLois Curfman McInnes smycols = mycols; 8908965ea79SLois Curfman McInnes svals = vals; 8918965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 8928965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 8938965ea79SLois Curfman McInnes smycols += ourlens[i]; 8948965ea79SLois Curfman McInnes svals += ourlens[i]; 8958965ea79SLois Curfman McInnes jj++; 8968965ea79SLois Curfman McInnes } 8978965ea79SLois Curfman McInnes 8988965ea79SLois Curfman McInnes /* read in other processors and ship out */ 8998965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 9008965ea79SLois Curfman McInnes nz = procsnz[i]; 9018965ea79SLois Curfman McInnes ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 9028965ea79SLois Curfman McInnes MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm); 9038965ea79SLois Curfman McInnes } 9048965ea79SLois Curfman McInnes PETSCFREE(procsnz); 9058965ea79SLois Curfman McInnes } 9068965ea79SLois Curfman McInnes else { 9078965ea79SLois Curfman McInnes /* receive numeric values */ 9088965ea79SLois Curfman McInnes vals = (Scalar*) PETSCMALLOC( nz*sizeof(Scalar) ); CHKPTRQ(vals); 9098965ea79SLois Curfman McInnes 9108965ea79SLois Curfman McInnes /* receive message of values*/ 9118965ea79SLois Curfman McInnes MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status); 9128965ea79SLois Curfman McInnes MPI_Get_count(&status,MPIU_SCALAR,&maxnz); 91339ddd567SLois Curfman McInnes if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file"); 9148965ea79SLois Curfman McInnes 9158965ea79SLois Curfman McInnes /* insert into matrix */ 9168965ea79SLois Curfman McInnes jj = rstart; 9178965ea79SLois Curfman McInnes smycols = mycols; 9188965ea79SLois Curfman McInnes svals = vals; 9198965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 9208965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 9218965ea79SLois Curfman McInnes smycols += ourlens[i]; 9228965ea79SLois Curfman McInnes svals += ourlens[i]; 9238965ea79SLois Curfman McInnes jj++; 9248965ea79SLois Curfman McInnes } 9258965ea79SLois Curfman McInnes } 9268965ea79SLois Curfman McInnes PETSCFREE(ourlens); PETSCFREE(vals); PETSCFREE(mycols); PETSCFREE(rowners); 9278965ea79SLois Curfman McInnes 9288965ea79SLois Curfman McInnes ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 9298965ea79SLois Curfman McInnes ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 9308965ea79SLois Curfman McInnes return 0; 9318965ea79SLois Curfman McInnes } 932