xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 4e220ebcc634eb0c8ccd78bed049803c9194c4d4)
18965ea79SLois Curfman McInnes #ifndef lint
2*4e220ebcSLois Curfman McInnes static char vcid[] = "$Id: mpidense.c,v 1.47 1996/08/15 16:05:25 bsmith Exp curfman $";
38965ea79SLois Curfman McInnes #endif
48965ea79SLois Curfman McInnes 
5ed3cc1f0SBarry Smith /*
6ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
7ed3cc1f0SBarry Smith */
8ed3cc1f0SBarry Smith 
970f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
10f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
118965ea79SLois Curfman McInnes 
127056b6fcSBarry Smith static int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
138965ea79SLois Curfman McInnes {
1439b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense *) mat->data;
1539b7565bSBarry Smith   int          ierr, i, j, rstart = A->rstart, rend = A->rend, row;
1639b7565bSBarry Smith   int          roworiented = A->roworiented;
178965ea79SLois Curfman McInnes 
1839b7565bSBarry Smith   if (A->insertmode != NOT_SET_VALUES && A->insertmode != addv) {
1939ddd567SLois Curfman McInnes     SETERRQ(1,"MatSetValues_MPIDense:Cannot mix inserts and adds");
208965ea79SLois Curfman McInnes   }
2139b7565bSBarry Smith   A->insertmode = addv;
228965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
2339ddd567SLois Curfman McInnes     if (idxm[i] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative row");
2439b7565bSBarry Smith     if (idxm[i] >= A->M) SETERRQ(1,"MatSetValues_MPIDense:Row too large");
258965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
268965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
2739b7565bSBarry Smith       if (roworiented) {
2839b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv); CHKERRQ(ierr);
2939b7565bSBarry Smith       }
3039b7565bSBarry Smith       else {
318965ea79SLois Curfman McInnes         for ( j=0; j<n; j++ ) {
3239ddd567SLois Curfman McInnes           if (idxn[j] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative column");
3339b7565bSBarry Smith           if (idxn[j] >= A->N) SETERRQ(1,"MatSetValues_MPIDense:Column too large");
3439b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv); CHKERRQ(ierr);
3539b7565bSBarry Smith         }
368965ea79SLois Curfman McInnes       }
378965ea79SLois Curfman McInnes     }
388965ea79SLois Curfman McInnes     else {
3939b7565bSBarry Smith       if (roworiented) {
4039b7565bSBarry Smith         ierr = StashValues_Private(&A->stash,idxm[i],n,idxn,v+i*n,addv); CHKERRQ(ierr);
4139b7565bSBarry Smith       }
4239b7565bSBarry Smith       else { /* must stash each seperately */
4339b7565bSBarry Smith         row = idxm[i];
4439b7565bSBarry Smith         for ( j=0; j<n; j++ ) {
457056b6fcSBarry Smith           ierr = StashValues_Private(&A->stash,row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
4639b7565bSBarry Smith         }
4739b7565bSBarry Smith       }
48b49de8d1SLois Curfman McInnes     }
49b49de8d1SLois Curfman McInnes   }
50b49de8d1SLois Curfman McInnes   return 0;
51b49de8d1SLois Curfman McInnes }
52b49de8d1SLois Curfman McInnes 
53b49de8d1SLois Curfman McInnes static int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
54b49de8d1SLois Curfman McInnes {
55b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
56b49de8d1SLois Curfman McInnes   int          ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row;
57b49de8d1SLois Curfman McInnes 
58b49de8d1SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
59b49de8d1SLois Curfman McInnes     if (idxm[i] < 0) SETERRQ(1,"MatGetValues_MPIDense:Negative row");
60b49de8d1SLois Curfman McInnes     if (idxm[i] >= mdn->M) SETERRQ(1,"MatGetValues_MPIDense:Row too large");
61b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
62b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
63b49de8d1SLois Curfman McInnes       for ( j=0; j<n; j++ ) {
64b49de8d1SLois Curfman McInnes         if (idxn[j] < 0) SETERRQ(1,"MatGetValues_MPIDense:Negative column");
65b49de8d1SLois Curfman McInnes         if (idxn[j] >= mdn->N)
66b49de8d1SLois Curfman McInnes           SETERRQ(1,"MatGetValues_MPIDense:Column too large");
67b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j); CHKERRQ(ierr);
68b49de8d1SLois Curfman McInnes       }
69b49de8d1SLois Curfman McInnes     }
70b49de8d1SLois Curfman McInnes     else {
71b49de8d1SLois Curfman McInnes       SETERRQ(1,"MatGetValues_MPIDense:Only local values currently supported");
728965ea79SLois Curfman McInnes     }
738965ea79SLois Curfman McInnes   }
748965ea79SLois Curfman McInnes   return 0;
758965ea79SLois Curfman McInnes }
768965ea79SLois Curfman McInnes 
77ff14e315SSatish Balay static int MatGetArray_MPIDense(Mat A,Scalar **array)
78ff14e315SSatish Balay {
79ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
80ff14e315SSatish Balay   int ierr;
81ff14e315SSatish Balay 
82ff14e315SSatish Balay   ierr = MatGetArray(a->A,array); CHKERRQ(ierr);
83ff14e315SSatish Balay   return 0;
84ff14e315SSatish Balay }
85ff14e315SSatish Balay 
86ff14e315SSatish Balay static int MatRestoreArray_MPIDense(Mat A,Scalar **array)
87ff14e315SSatish Balay {
88ff14e315SSatish Balay   return 0;
89ff14e315SSatish Balay }
90ff14e315SSatish Balay 
9139ddd567SLois Curfman McInnes static int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
928965ea79SLois Curfman McInnes {
9339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
948965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
9539ddd567SLois Curfman McInnes   int          size = mdn->size, *owners = mdn->rowners, rank = mdn->rank;
968965ea79SLois Curfman McInnes   int          *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
9739ddd567SLois Curfman McInnes   int          tag = mat->tag, *owner,*starts,count,ierr;
988965ea79SLois Curfman McInnes   InsertMode   addv;
9939ddd567SLois Curfman McInnes   MPI_Request  *send_waits,*recv_waits;
1008965ea79SLois Curfman McInnes   Scalar       *rvalues,*svalues;
1018965ea79SLois Curfman McInnes 
1028965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
1036d84be18SBarry Smith   MPI_Allreduce(&mdn->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
1047056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
1057056b6fcSBarry Smith     SETERRQ(1,"MatAssemblyBegin_MPIDense:Cannot mix adds/inserts on different procs");
1068965ea79SLois Curfman McInnes   }
10739ddd567SLois Curfman McInnes   mdn->insertmode = addv; /* in case this processor had no cache */
1088965ea79SLois Curfman McInnes 
1098965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
1100452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
111cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
1120452661fSBarry Smith   owner = (int *) PetscMalloc( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
11339ddd567SLois Curfman McInnes   for ( i=0; i<mdn->stash.n; i++ ) {
11439ddd567SLois Curfman McInnes     idx = mdn->stash.idx[i];
1158965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
1168965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
1178965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
1188965ea79SLois Curfman McInnes       }
1198965ea79SLois Curfman McInnes     }
1208965ea79SLois Curfman McInnes   }
1218965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
1228965ea79SLois Curfman McInnes 
1238965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
1240452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
1256d84be18SBarry Smith   MPI_Allreduce(procs,work,size,MPI_INT,MPI_SUM,comm);
1268965ea79SLois Curfman McInnes   nreceives = work[rank];
1276d84be18SBarry Smith   MPI_Allreduce(nprocs,work,size,MPI_INT,MPI_MAX,comm);
1288965ea79SLois Curfman McInnes   nmax = work[rank];
1290452661fSBarry Smith   PetscFree(work);
1308965ea79SLois Curfman McInnes 
1318965ea79SLois Curfman McInnes   /* post receives:
1328965ea79SLois Curfman McInnes        1) each message will consist of ordered pairs
1338965ea79SLois Curfman McInnes      (global index,value) we store the global index as a double
1348965ea79SLois Curfman McInnes      to simplify the message passing.
1358965ea79SLois Curfman McInnes        2) since we don't know how long each individual message is we
1368965ea79SLois Curfman McInnes      allocate the largest needed buffer for each receive. Potentially
1378965ea79SLois Curfman McInnes      this is a lot of wasted space.
1388965ea79SLois Curfman McInnes 
1398965ea79SLois Curfman McInnes        This could be done better.
1408965ea79SLois Curfman McInnes   */
1410452661fSBarry Smith   rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));
1428965ea79SLois Curfman McInnes   CHKPTRQ(rvalues);
1430452661fSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));
1448965ea79SLois Curfman McInnes   CHKPTRQ(recv_waits);
1458965ea79SLois Curfman McInnes   for ( i=0; i<nreceives; i++ ) {
1467056b6fcSBarry Smith     MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
1478965ea79SLois Curfman McInnes   }
1488965ea79SLois Curfman McInnes 
1498965ea79SLois Curfman McInnes   /* do sends:
1508965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
1518965ea79SLois Curfman McInnes          the ith processor
1528965ea79SLois Curfman McInnes   */
1530452661fSBarry Smith   svalues = (Scalar *) PetscMalloc( 3*(mdn->stash.n+1)*sizeof(Scalar) );
15439ddd567SLois Curfman McInnes   CHKPTRQ(svalues);
1550452661fSBarry Smith   send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));
1568965ea79SLois Curfman McInnes   CHKPTRQ(send_waits);
1570452661fSBarry Smith   starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts);
1588965ea79SLois Curfman McInnes   starts[0] = 0;
1598965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
16039ddd567SLois Curfman McInnes   for ( i=0; i<mdn->stash.n; i++ ) {
16139ddd567SLois Curfman McInnes     svalues[3*starts[owner[i]]]       = (Scalar)  mdn->stash.idx[i];
16239ddd567SLois Curfman McInnes     svalues[3*starts[owner[i]]+1]     = (Scalar)  mdn->stash.idy[i];
16339ddd567SLois Curfman McInnes     svalues[3*(starts[owner[i]]++)+2] =  mdn->stash.array[i];
1648965ea79SLois Curfman McInnes   }
1650452661fSBarry Smith   PetscFree(owner);
1668965ea79SLois Curfman McInnes   starts[0] = 0;
1678965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1688965ea79SLois Curfman McInnes   count = 0;
1698965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
1708965ea79SLois Curfman McInnes     if (procs[i]) {
1717056b6fcSBarry Smith       MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);
1728965ea79SLois Curfman McInnes     }
1738965ea79SLois Curfman McInnes   }
1740452661fSBarry Smith   PetscFree(starts); PetscFree(nprocs);
1758965ea79SLois Curfman McInnes 
1768965ea79SLois Curfman McInnes   /* Free cache space */
17739ddd567SLois Curfman McInnes   ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr);
1788965ea79SLois Curfman McInnes 
17939ddd567SLois Curfman McInnes   mdn->svalues    = svalues;    mdn->rvalues = rvalues;
18039ddd567SLois Curfman McInnes   mdn->nsends     = nsends;     mdn->nrecvs = nreceives;
18139ddd567SLois Curfman McInnes   mdn->send_waits = send_waits; mdn->recv_waits = recv_waits;
18239ddd567SLois Curfman McInnes   mdn->rmax       = nmax;
1838965ea79SLois Curfman McInnes 
1848965ea79SLois Curfman McInnes   return 0;
1858965ea79SLois Curfman McInnes }
18639ddd567SLois Curfman McInnes extern int MatSetUpMultiply_MPIDense(Mat);
1878965ea79SLois Curfman McInnes 
18839ddd567SLois Curfman McInnes static int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
1898965ea79SLois Curfman McInnes {
19039ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
1918965ea79SLois Curfman McInnes   MPI_Status   *send_status,recv_status;
19239ddd567SLois Curfman McInnes   int          imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col;
1938965ea79SLois Curfman McInnes   Scalar       *values,val;
19439ddd567SLois Curfman McInnes   InsertMode   addv = mdn->insertmode;
1958965ea79SLois Curfman McInnes 
1968965ea79SLois Curfman McInnes   /*  wait on receives */
1978965ea79SLois Curfman McInnes   while (count) {
19839ddd567SLois Curfman McInnes     MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status);
1998965ea79SLois Curfman McInnes     /* unpack receives into our local space */
20039ddd567SLois Curfman McInnes     values = mdn->rvalues + 3*imdex*mdn->rmax;
2018965ea79SLois Curfman McInnes     MPI_Get_count(&recv_status,MPIU_SCALAR,&n);
2028965ea79SLois Curfman McInnes     n = n/3;
2038965ea79SLois Curfman McInnes     for ( i=0; i<n; i++ ) {
204227d817aSBarry Smith       row = (int) PetscReal(values[3*i]) - mdn->rstart;
205227d817aSBarry Smith       col = (int) PetscReal(values[3*i+1]);
2068965ea79SLois Curfman McInnes       val = values[3*i+2];
20739ddd567SLois Curfman McInnes       if (col >= 0 && col < mdn->N) {
20839ddd567SLois Curfman McInnes         MatSetValues(mdn->A,1,&row,1,&col,&val,addv);
2098965ea79SLois Curfman McInnes       }
21039ddd567SLois Curfman McInnes       else {SETERRQ(1,"MatAssemblyEnd_MPIDense:Invalid column");}
2118965ea79SLois Curfman McInnes     }
2128965ea79SLois Curfman McInnes     count--;
2138965ea79SLois Curfman McInnes   }
2140452661fSBarry Smith   PetscFree(mdn->recv_waits); PetscFree(mdn->rvalues);
2158965ea79SLois Curfman McInnes 
2168965ea79SLois Curfman McInnes   /* wait on sends */
21739ddd567SLois Curfman McInnes   if (mdn->nsends) {
2187056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(mdn->nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
21939ddd567SLois Curfman McInnes     MPI_Waitall(mdn->nsends,mdn->send_waits,send_status);
2200452661fSBarry Smith     PetscFree(send_status);
2218965ea79SLois Curfman McInnes   }
2220452661fSBarry Smith   PetscFree(mdn->send_waits); PetscFree(mdn->svalues);
2238965ea79SLois Curfman McInnes 
22439ddd567SLois Curfman McInnes   mdn->insertmode = NOT_SET_VALUES;
22539ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr);
22639ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr);
2278965ea79SLois Curfman McInnes 
2286d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
22939ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr);
2308965ea79SLois Curfman McInnes   }
2318965ea79SLois Curfman McInnes   return 0;
2328965ea79SLois Curfman McInnes }
2338965ea79SLois Curfman McInnes 
23439ddd567SLois Curfman McInnes static int MatZeroEntries_MPIDense(Mat A)
2358965ea79SLois Curfman McInnes {
23639ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense *) A->data;
23739ddd567SLois Curfman McInnes   return MatZeroEntries(l->A);
2388965ea79SLois Curfman McInnes }
2398965ea79SLois Curfman McInnes 
240*4e220ebcSLois Curfman McInnes static int MatGetBlockSize_MPIDense(Mat A,int *bs)
241*4e220ebcSLois Curfman McInnes {
242*4e220ebcSLois Curfman McInnes   *bs = 1;
243*4e220ebcSLois Curfman McInnes   return 0;
244*4e220ebcSLois Curfman McInnes }
245*4e220ebcSLois Curfman McInnes 
2468965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2478965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2488965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2493501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2508965ea79SLois Curfman McInnes    routine.
2518965ea79SLois Curfman McInnes */
25239ddd567SLois Curfman McInnes static int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag)
2538965ea79SLois Curfman McInnes {
25439ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense *) A->data;
2558965ea79SLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
2568965ea79SLois Curfman McInnes   int            *procs,*nprocs,j,found,idx,nsends,*work;
2578965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2588965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2598965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2608965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2618965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2628965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2638965ea79SLois Curfman McInnes   IS             istmp;
2648965ea79SLois Curfman McInnes 
26577c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
2668965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
2678965ea79SLois Curfman McInnes 
2688965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
2690452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
270cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
2710452661fSBarry Smith   owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
2728965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
2738965ea79SLois Curfman McInnes     idx = rows[i];
2748965ea79SLois Curfman McInnes     found = 0;
2758965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
2768965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
2778965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
2788965ea79SLois Curfman McInnes       }
2798965ea79SLois Curfman McInnes     }
28039ddd567SLois Curfman McInnes     if (!found) SETERRQ(1,"MatZeroRows_MPIDense:Index out of range");
2818965ea79SLois Curfman McInnes   }
2828965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
2838965ea79SLois Curfman McInnes 
2848965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
2850452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
2868965ea79SLois Curfman McInnes   MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);
2878965ea79SLois Curfman McInnes   nrecvs = work[rank];
2888965ea79SLois Curfman McInnes   MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);
2898965ea79SLois Curfman McInnes   nmax = work[rank];
2900452661fSBarry Smith   PetscFree(work);
2918965ea79SLois Curfman McInnes 
2928965ea79SLois Curfman McInnes   /* post receives:   */
2930452661fSBarry Smith   rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */
2948965ea79SLois Curfman McInnes   CHKPTRQ(rvalues);
2950452661fSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));
2968965ea79SLois Curfman McInnes   CHKPTRQ(recv_waits);
2978965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
2988965ea79SLois Curfman McInnes     MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
2998965ea79SLois Curfman McInnes   }
3008965ea79SLois Curfman McInnes 
3018965ea79SLois Curfman McInnes   /* do sends:
3028965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3038965ea79SLois Curfman McInnes          the ith processor
3048965ea79SLois Curfman McInnes   */
3050452661fSBarry Smith   svalues    = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
3067056b6fcSBarry Smith   send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
3070452661fSBarry Smith   starts     = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts);
3088965ea79SLois Curfman McInnes   starts[0]  = 0;
3098965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3108965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
3118965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3128965ea79SLois Curfman McInnes   }
3138965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3148965ea79SLois Curfman McInnes 
3158965ea79SLois Curfman McInnes   starts[0] = 0;
3168965ea79SLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3178965ea79SLois Curfman McInnes   count = 0;
3188965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3198965ea79SLois Curfman McInnes     if (procs[i]) {
3208965ea79SLois Curfman McInnes       MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);
3218965ea79SLois Curfman McInnes     }
3228965ea79SLois Curfman McInnes   }
3230452661fSBarry Smith   PetscFree(starts);
3248965ea79SLois Curfman McInnes 
3258965ea79SLois Curfman McInnes   base = owners[rank];
3268965ea79SLois Curfman McInnes 
3278965ea79SLois Curfman McInnes   /*  wait on receives */
3280452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
3298965ea79SLois Curfman McInnes   source = lens + nrecvs;
3308965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3318965ea79SLois Curfman McInnes   while (count) {
3328965ea79SLois Curfman McInnes     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
3338965ea79SLois Curfman McInnes     /* unpack receives into our local space */
3348965ea79SLois Curfman McInnes     MPI_Get_count(&recv_status,MPI_INT,&n);
3358965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3368965ea79SLois Curfman McInnes     lens[imdex]  = n;
3378965ea79SLois Curfman McInnes     slen += n;
3388965ea79SLois Curfman McInnes     count--;
3398965ea79SLois Curfman McInnes   }
3400452661fSBarry Smith   PetscFree(recv_waits);
3418965ea79SLois Curfman McInnes 
3428965ea79SLois Curfman McInnes   /* move the data into the send scatter */
3430452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
3448965ea79SLois Curfman McInnes   count = 0;
3458965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
3468965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3478965ea79SLois Curfman McInnes     for ( j=0; j<lens[i]; j++ ) {
3488965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3498965ea79SLois Curfman McInnes     }
3508965ea79SLois Curfman McInnes   }
3510452661fSBarry Smith   PetscFree(rvalues); PetscFree(lens);
3520452661fSBarry Smith   PetscFree(owner); PetscFree(nprocs);
3538965ea79SLois Curfman McInnes 
3548965ea79SLois Curfman McInnes   /* actually zap the local rows */
355537820f0SBarry Smith   ierr = ISCreateGeneral(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
3568965ea79SLois Curfman McInnes   PLogObjectParent(A,istmp);
3570452661fSBarry Smith   PetscFree(lrows);
3588965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
3598965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp); CHKERRQ(ierr);
3608965ea79SLois Curfman McInnes 
3618965ea79SLois Curfman McInnes   /* wait on sends */
3628965ea79SLois Curfman McInnes   if (nsends) {
3637056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
3648965ea79SLois Curfman McInnes     MPI_Waitall(nsends,send_waits,send_status);
3650452661fSBarry Smith     PetscFree(send_status);
3668965ea79SLois Curfman McInnes   }
3670452661fSBarry Smith   PetscFree(send_waits); PetscFree(svalues);
3688965ea79SLois Curfman McInnes 
3698965ea79SLois Curfman McInnes   return 0;
3708965ea79SLois Curfman McInnes }
3718965ea79SLois Curfman McInnes 
37239ddd567SLois Curfman McInnes static int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3738965ea79SLois Curfman McInnes {
37439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
3758965ea79SLois Curfman McInnes   int          ierr;
376c456f294SBarry Smith 
3777056b6fcSBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
3787056b6fcSBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
37944cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy); CHKERRQ(ierr);
3808965ea79SLois Curfman McInnes   return 0;
3818965ea79SLois Curfman McInnes }
3828965ea79SLois Curfman McInnes 
38339ddd567SLois Curfman McInnes static int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
3848965ea79SLois Curfman McInnes {
38539ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
3868965ea79SLois Curfman McInnes   int          ierr;
387c456f294SBarry Smith 
388c456f294SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
389c456f294SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
39044cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr);
3918965ea79SLois Curfman McInnes   return 0;
3928965ea79SLois Curfman McInnes }
3938965ea79SLois Curfman McInnes 
394096963f5SLois Curfman McInnes static int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy)
395096963f5SLois Curfman McInnes {
396096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
397096963f5SLois Curfman McInnes   int          ierr;
3983501a2bdSLois Curfman McInnes   Scalar       zero = 0.0;
399096963f5SLois Curfman McInnes 
4003501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy); CHKERRQ(ierr);
40144cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr);
402537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
403537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
404096963f5SLois Curfman McInnes   return 0;
405096963f5SLois Curfman McInnes }
406096963f5SLois Curfman McInnes 
407096963f5SLois Curfman McInnes static int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
408096963f5SLois Curfman McInnes {
409096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
410096963f5SLois Curfman McInnes   int          ierr;
411096963f5SLois Curfman McInnes 
4123501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz); CHKERRQ(ierr);
41344cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr);
414537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
415537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
416096963f5SLois Curfman McInnes   return 0;
417096963f5SLois Curfman McInnes }
418096963f5SLois Curfman McInnes 
41939ddd567SLois Curfman McInnes static int MatGetDiagonal_MPIDense(Mat A,Vec v)
4208965ea79SLois Curfman McInnes {
42139ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
422096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data;
42344cd7ae7SLois Curfman McInnes   int          ierr, len, i, n, m = a->m, radd;
42444cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
425ed3cc1f0SBarry Smith 
42644cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
427096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x); CHKERRQ(ierr);
428096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n); CHKERRQ(ierr);
429096963f5SLois Curfman McInnes   if (n != a->M) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
43044cd7ae7SLois Curfman McInnes   len = PetscMin(aloc->m,aloc->n);
4317ddc982cSLois Curfman McInnes   radd = a->rstart*m;
43244cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
433096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
434096963f5SLois Curfman McInnes   }
435096963f5SLois Curfman McInnes   return 0;
4368965ea79SLois Curfman McInnes }
4378965ea79SLois Curfman McInnes 
43839ddd567SLois Curfman McInnes static int MatDestroy_MPIDense(PetscObject obj)
4398965ea79SLois Curfman McInnes {
4408965ea79SLois Curfman McInnes   Mat          mat = (Mat) obj;
4413501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4428965ea79SLois Curfman McInnes   int          ierr;
443ed3cc1f0SBarry Smith 
4448965ea79SLois Curfman McInnes #if defined(PETSC_LOG)
4453501a2bdSLois Curfman McInnes   PLogObjectState(obj,"Rows=%d, Cols=%d",mdn->M,mdn->N);
4468965ea79SLois Curfman McInnes #endif
4470452661fSBarry Smith   PetscFree(mdn->rowners);
4483501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A); CHKERRQ(ierr);
4493501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4503501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
451622d7880SLois Curfman McInnes   if (mdn->factor) {
452622d7880SLois Curfman McInnes     if (mdn->factor->temp)   PetscFree(mdn->factor->temp);
453622d7880SLois Curfman McInnes     if (mdn->factor->tag)    PetscFree(mdn->factor->tag);
454622d7880SLois Curfman McInnes     if (mdn->factor->pivots) PetscFree(mdn->factor->pivots);
455622d7880SLois Curfman McInnes     PetscFree(mdn->factor);
456622d7880SLois Curfman McInnes   }
4570452661fSBarry Smith   PetscFree(mdn);
4588965ea79SLois Curfman McInnes   PLogObjectDestroy(mat);
4590452661fSBarry Smith   PetscHeaderDestroy(mat);
4608965ea79SLois Curfman McInnes   return 0;
4618965ea79SLois Curfman McInnes }
46239ddd567SLois Curfman McInnes 
4638965ea79SLois Curfman McInnes #include "pinclude/pviewer.h"
4648965ea79SLois Curfman McInnes 
46539ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
4668965ea79SLois Curfman McInnes {
46739ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4688965ea79SLois Curfman McInnes   int          ierr;
4697056b6fcSBarry Smith 
47039ddd567SLois Curfman McInnes   if (mdn->size == 1) {
47139ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
4728965ea79SLois Curfman McInnes   }
47339ddd567SLois Curfman McInnes   else SETERRQ(1,"MatView_MPIDense_Binary:Only uniprocessor output supported");
4748965ea79SLois Curfman McInnes   return 0;
4758965ea79SLois Curfman McInnes }
4768965ea79SLois Curfman McInnes 
47739ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer)
4788965ea79SLois Curfman McInnes {
47939ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
48039ddd567SLois Curfman McInnes   int          ierr, format;
4818965ea79SLois Curfman McInnes   FILE         *fd;
48219bcc07fSBarry Smith   ViewerType   vtype;
4838965ea79SLois Curfman McInnes 
48419bcc07fSBarry Smith   ViewerGetType(viewer,&vtype);
48590ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
48690ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
48790ace30eSBarry Smith   if (format == ASCII_FORMAT_INFO_DETAILED) {
488*4e220ebcSLois Curfman McInnes     int rank;
489*4e220ebcSLois Curfman McInnes     MatInfo info;
490096963f5SLois Curfman McInnes     MPI_Comm_rank(mat->comm,&rank);
491*4e220ebcSLois Curfman McInnes     ierr = MatGetInfo(mat,MAT_LOCAL,&info);
49277c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
493*4e220ebcSLois Curfman McInnes       fprintf(fd,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m,
494*4e220ebcSLois Curfman McInnes          (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
495096963f5SLois Curfman McInnes       fflush(fd);
49677c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
4973501a2bdSLois Curfman McInnes     ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr);
4983501a2bdSLois Curfman McInnes     return 0;
4993501a2bdSLois Curfman McInnes   }
50090ace30eSBarry Smith   else if (format == ASCII_FORMAT_INFO) {
501096963f5SLois Curfman McInnes     return 0;
5028965ea79SLois Curfman McInnes   }
50319bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
50477c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
50539ddd567SLois Curfman McInnes     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n",
50639ddd567SLois Curfman McInnes              mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n);
50739ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5088965ea79SLois Curfman McInnes     fflush(fd);
50977c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
5108965ea79SLois Curfman McInnes   }
5118965ea79SLois Curfman McInnes   else {
51239ddd567SLois Curfman McInnes     int size = mdn->size, rank = mdn->rank;
5138965ea79SLois Curfman McInnes     if (size == 1) {
51439ddd567SLois Curfman McInnes       ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5158965ea79SLois Curfman McInnes     }
5168965ea79SLois Curfman McInnes     else {
5178965ea79SLois Curfman McInnes       /* assemble the entire matrix onto first processor. */
5188965ea79SLois Curfman McInnes       Mat          A;
51939ddd567SLois Curfman McInnes       int          M = mdn->M, N = mdn->N,m,row,i, nz, *cols;
52039ddd567SLois Curfman McInnes       Scalar       *vals;
52139ddd567SLois Curfman McInnes       Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data;
5228965ea79SLois Curfman McInnes 
5238965ea79SLois Curfman McInnes       if (!rank) {
524b4fd4287SBarry Smith         ierr = MatCreateMPIDense(mat->comm,M,M,N,N,PETSC_NULL,&A); CHKERRQ(ierr);
5258965ea79SLois Curfman McInnes       }
5268965ea79SLois Curfman McInnes       else {
527b4fd4287SBarry Smith         ierr = MatCreateMPIDense(mat->comm,0,M,N,N,PETSC_NULL,&A); CHKERRQ(ierr);
5288965ea79SLois Curfman McInnes       }
5298965ea79SLois Curfman McInnes       PLogObjectParent(mat,A);
5308965ea79SLois Curfman McInnes 
53139ddd567SLois Curfman McInnes       /* Copy the matrix ... This isn't the most efficient means,
53239ddd567SLois Curfman McInnes          but it's quick for now */
53339ddd567SLois Curfman McInnes       row = mdn->rstart; m = Amdn->m;
5348965ea79SLois Curfman McInnes       for ( i=0; i<m; i++ ) {
53539ddd567SLois Curfman McInnes         ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
53639ddd567SLois Curfman McInnes         ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr);
53739ddd567SLois Curfman McInnes         ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
53839ddd567SLois Curfman McInnes         row++;
5398965ea79SLois Curfman McInnes       }
5408965ea79SLois Curfman McInnes 
5416d4a8577SBarry Smith       ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5426d4a8577SBarry Smith       ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5438965ea79SLois Curfman McInnes       if (!rank) {
54439ddd567SLois Curfman McInnes         ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr);
5458965ea79SLois Curfman McInnes       }
5468965ea79SLois Curfman McInnes       ierr = MatDestroy(A); CHKERRQ(ierr);
5478965ea79SLois Curfman McInnes     }
5488965ea79SLois Curfman McInnes   }
5498965ea79SLois Curfman McInnes   return 0;
5508965ea79SLois Curfman McInnes }
5518965ea79SLois Curfman McInnes 
55239ddd567SLois Curfman McInnes static int MatView_MPIDense(PetscObject obj,Viewer viewer)
5538965ea79SLois Curfman McInnes {
5548965ea79SLois Curfman McInnes   Mat          mat = (Mat) obj;
55539ddd567SLois Curfman McInnes   int          ierr;
556bcd2baecSBarry Smith   ViewerType   vtype;
5578965ea79SLois Curfman McInnes 
558bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
559bcd2baecSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
56039ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
5618965ea79SLois Curfman McInnes   }
562bcd2baecSBarry Smith   else if (vtype == ASCII_FILES_VIEWER) {
56339ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
5648965ea79SLois Curfman McInnes   }
565bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
56639ddd567SLois Curfman McInnes     return MatView_MPIDense_Binary(mat,viewer);
5678965ea79SLois Curfman McInnes   }
5688965ea79SLois Curfman McInnes   return 0;
5698965ea79SLois Curfman McInnes }
5708965ea79SLois Curfman McInnes 
571*4e220ebcSLois Curfman McInnes static int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
5728965ea79SLois Curfman McInnes {
5733501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
5743501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
575*4e220ebcSLois Curfman McInnes   int          ierr;
576*4e220ebcSLois Curfman McInnes   double       isend[5], irecv[5];
5778965ea79SLois Curfman McInnes 
578*4e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
579*4e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
580*4e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
581*4e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
582*4e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
583*4e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info); CHKERRQ(ierr);
584*4e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
585*4e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
5868965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
587*4e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
588*4e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
589*4e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
590*4e220ebcSLois Curfman McInnes     info->memory       = isend[3];
591*4e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
5928965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
5933501a2bdSLois Curfman McInnes     MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_MAX,A->comm);
594*4e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
595*4e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
596*4e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
597*4e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
598*4e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
5998965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
6003501a2bdSLois Curfman McInnes     MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_SUM,A->comm);
601*4e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
602*4e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
603*4e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
604*4e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
605*4e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6068965ea79SLois Curfman McInnes   }
607*4e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
608*4e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
609*4e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6108965ea79SLois Curfman McInnes   return 0;
6118965ea79SLois Curfman McInnes }
6128965ea79SLois Curfman McInnes 
6138c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*);
6148aaee692SLois Curfman McInnes    extern int MatLUFactorNumeric_MPIDense(Mat,Mat*);
6158aaee692SLois Curfman McInnes    extern int MatLUFactor_MPIDense(Mat,IS,IS,double);
6168aaee692SLois Curfman McInnes    extern int MatSolve_MPIDense(Mat,Vec,Vec);
6178c469469SLois Curfman McInnes    extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec);
6188aaee692SLois Curfman McInnes    extern int MatSolveTrans_MPIDense(Mat,Vec,Vec);
6198aaee692SLois Curfman McInnes    extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */
6208aaee692SLois Curfman McInnes 
62139ddd567SLois Curfman McInnes static int MatSetOption_MPIDense(Mat A,MatOption op)
6228965ea79SLois Curfman McInnes {
62339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
6248965ea79SLois Curfman McInnes 
6256d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6266d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
6276d4a8577SBarry Smith       op == MAT_COLUMNS_SORTED ||
6286d4a8577SBarry Smith       op == MAT_ROW_ORIENTED) {
6298965ea79SLois Curfman McInnes         MatSetOption(a->A,op);
6308965ea79SLois Curfman McInnes   }
6316d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
6326d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
6336d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
6346d4a8577SBarry Smith            op == MAT_YES_NEW_DIAGONALS)
63594a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_MPIDense:Option ignored\n");
6366d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)
63739b7565bSBarry Smith     { a->roworiented = 0; MatSetOption(a->A,op);}
6386d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
6396d4a8577SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:MAT_NO_NEW_DIAGONALS");}
6408965ea79SLois Curfman McInnes   else
64139ddd567SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:unknown option");}
6428965ea79SLois Curfman McInnes   return 0;
6438965ea79SLois Curfman McInnes }
6448965ea79SLois Curfman McInnes 
6453501a2bdSLois Curfman McInnes static int MatGetSize_MPIDense(Mat A,int *m,int *n)
6468965ea79SLois Curfman McInnes {
6473501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6488965ea79SLois Curfman McInnes   *m = mat->M; *n = mat->N;
6498965ea79SLois Curfman McInnes   return 0;
6508965ea79SLois Curfman McInnes }
6518965ea79SLois Curfman McInnes 
6523501a2bdSLois Curfman McInnes static int MatGetLocalSize_MPIDense(Mat A,int *m,int *n)
6538965ea79SLois Curfman McInnes {
6543501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6558965ea79SLois Curfman McInnes   *m = mat->m; *n = mat->N;
6568965ea79SLois Curfman McInnes   return 0;
6578965ea79SLois Curfman McInnes }
6588965ea79SLois Curfman McInnes 
6593501a2bdSLois Curfman McInnes static int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
6608965ea79SLois Curfman McInnes {
6613501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6628965ea79SLois Curfman McInnes   *m = mat->rstart; *n = mat->rend;
6638965ea79SLois Curfman McInnes   return 0;
6648965ea79SLois Curfman McInnes }
6658965ea79SLois Curfman McInnes 
6663501a2bdSLois Curfman McInnes static int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
6678965ea79SLois Curfman McInnes {
6683501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
66939ddd567SLois Curfman McInnes   int          lrow, rstart = mat->rstart, rend = mat->rend;
6708965ea79SLois Curfman McInnes 
67139ddd567SLois Curfman McInnes   if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIDense:only local rows")
6728965ea79SLois Curfman McInnes   lrow = row - rstart;
67339ddd567SLois Curfman McInnes   return MatGetRow(mat->A,lrow,nz,idx,v);
6748965ea79SLois Curfman McInnes }
6758965ea79SLois Curfman McInnes 
67639ddd567SLois Curfman McInnes static int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
6778965ea79SLois Curfman McInnes {
6780452661fSBarry Smith   if (idx) PetscFree(*idx);
6790452661fSBarry Smith   if (v) PetscFree(*v);
6808965ea79SLois Curfman McInnes   return 0;
6818965ea79SLois Curfman McInnes }
6828965ea79SLois Curfman McInnes 
683cddf8d76SBarry Smith static int MatNorm_MPIDense(Mat A,NormType type,double *norm)
684096963f5SLois Curfman McInnes {
6853501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
6863501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
6873501a2bdSLois Curfman McInnes   int          ierr, i, j;
6883501a2bdSLois Curfman McInnes   double       sum = 0.0;
6893501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
6903501a2bdSLois Curfman McInnes 
6913501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
6923501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm); CHKERRQ(ierr);
6933501a2bdSLois Curfman McInnes   } else {
6943501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
6953501a2bdSLois Curfman McInnes       for (i=0; i<mat->n*mat->m; i++ ) {
6963501a2bdSLois Curfman McInnes #if defined(PETSC_COMPLEX)
6973501a2bdSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
6983501a2bdSLois Curfman McInnes #else
6993501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7003501a2bdSLois Curfman McInnes #endif
7013501a2bdSLois Curfman McInnes       }
7026d84be18SBarry Smith       MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);
7033501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
7043501a2bdSLois Curfman McInnes       PLogFlops(2*mat->n*mat->m);
7053501a2bdSLois Curfman McInnes     }
7063501a2bdSLois Curfman McInnes     else if (type == NORM_1) {
7073501a2bdSLois Curfman McInnes       double *tmp, *tmp2;
7080452661fSBarry Smith       tmp  = (double *) PetscMalloc( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp);
7093501a2bdSLois Curfman McInnes       tmp2 = tmp + mdn->N;
710cddf8d76SBarry Smith       PetscMemzero(tmp,2*mdn->N*sizeof(double));
711096963f5SLois Curfman McInnes       *norm = 0.0;
7123501a2bdSLois Curfman McInnes       v = mat->v;
7133501a2bdSLois Curfman McInnes       for ( j=0; j<mat->n; j++ ) {
7143501a2bdSLois Curfman McInnes         for ( i=0; i<mat->m; i++ ) {
71567e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
7163501a2bdSLois Curfman McInnes         }
7173501a2bdSLois Curfman McInnes       }
7186d84be18SBarry Smith       MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);
7193501a2bdSLois Curfman McInnes       for ( j=0; j<mdn->N; j++ ) {
7203501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
7213501a2bdSLois Curfman McInnes       }
7220452661fSBarry Smith       PetscFree(tmp);
7233501a2bdSLois Curfman McInnes       PLogFlops(mat->n*mat->m);
7243501a2bdSLois Curfman McInnes     }
7253501a2bdSLois Curfman McInnes     else if (type == NORM_INFINITY) { /* max row norm */
7263501a2bdSLois Curfman McInnes       double ntemp;
7273501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr);
7286d84be18SBarry Smith       MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);
7293501a2bdSLois Curfman McInnes     }
7303501a2bdSLois Curfman McInnes     else {
7313501a2bdSLois Curfman McInnes       SETERRQ(1,"MatNorm_MPIDense:No support for two norm");
7323501a2bdSLois Curfman McInnes     }
7333501a2bdSLois Curfman McInnes   }
7343501a2bdSLois Curfman McInnes   return 0;
7353501a2bdSLois Curfman McInnes }
7363501a2bdSLois Curfman McInnes 
7373501a2bdSLois Curfman McInnes static int MatTranspose_MPIDense(Mat A,Mat *matout)
7383501a2bdSLois Curfman McInnes {
7393501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
7403501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data;
7413501a2bdSLois Curfman McInnes   Mat          B;
7423501a2bdSLois Curfman McInnes   int          M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart;
7433501a2bdSLois Curfman McInnes   int          j, i, ierr;
7443501a2bdSLois Curfman McInnes   Scalar       *v;
7453501a2bdSLois Curfman McInnes 
7467056b6fcSBarry Smith   if (matout == PETSC_NULL && M != N) {
7473501a2bdSLois Curfman McInnes     SETERRQ(1,"MatTranspose_MPIDense:Supports square matrix only in-place");
7487056b6fcSBarry Smith   }
7497056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
7503501a2bdSLois Curfman McInnes 
7513501a2bdSLois Curfman McInnes   m = Aloc->m; n = Aloc->n; v = Aloc->v;
7520452661fSBarry Smith   rwork = (int *) PetscMalloc(n*sizeof(int)); CHKPTRQ(rwork);
7533501a2bdSLois Curfman McInnes   for ( j=0; j<n; j++ ) {
7543501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
7553501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr);
7563501a2bdSLois Curfman McInnes     v   += m;
7573501a2bdSLois Curfman McInnes   }
7580452661fSBarry Smith   PetscFree(rwork);
7596d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
7606d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
7613638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
7623501a2bdSLois Curfman McInnes     *matout = B;
7633501a2bdSLois Curfman McInnes   } else {
7643501a2bdSLois Curfman McInnes     /* This isn't really an in-place transpose, but free data struct from a */
7650452661fSBarry Smith     PetscFree(a->rowners);
7663501a2bdSLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
7673501a2bdSLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
7683501a2bdSLois Curfman McInnes     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
7690452661fSBarry Smith     PetscFree(a);
7703501a2bdSLois Curfman McInnes     PetscMemcpy(A,B,sizeof(struct _Mat));
7710452661fSBarry Smith     PetscHeaderDestroy(B);
7723501a2bdSLois Curfman McInnes   }
773096963f5SLois Curfman McInnes   return 0;
774096963f5SLois Curfman McInnes }
775096963f5SLois Curfman McInnes 
77644cd7ae7SLois Curfman McInnes #include "pinclude/plapack.h"
77744cd7ae7SLois Curfman McInnes static int MatScale_MPIDense(Scalar *alpha,Mat inA)
77844cd7ae7SLois Curfman McInnes {
77944cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense *) inA->data;
78044cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->A->data;
78144cd7ae7SLois Curfman McInnes   int          one = 1, nz;
78244cd7ae7SLois Curfman McInnes 
78344cd7ae7SLois Curfman McInnes   nz = a->m*a->n;
78444cd7ae7SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
78544cd7ae7SLois Curfman McInnes   PLogFlops(nz);
78644cd7ae7SLois Curfman McInnes   return 0;
78744cd7ae7SLois Curfman McInnes }
78844cd7ae7SLois Curfman McInnes 
7893d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat,Mat *,int);
7908965ea79SLois Curfman McInnes 
7918965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
79239ddd567SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIDense,
79339ddd567SLois Curfman McInnes        MatGetRow_MPIDense,MatRestoreRow_MPIDense,
79439ddd567SLois Curfman McInnes        MatMult_MPIDense,MatMultAdd_MPIDense,
795096963f5SLois Curfman McInnes        MatMultTrans_MPIDense,MatMultTransAdd_MPIDense,
796e7ca0642SLois Curfman McInnes /*       MatSolve_MPIDense,0, */
797e7ca0642SLois Curfman McInnes        0,0,
79839ddd567SLois Curfman McInnes        0,0,
7998c469469SLois Curfman McInnes        0,0,
8008c469469SLois Curfman McInnes /*       MatLUFactor_MPIDense,0, */
8018aaee692SLois Curfman McInnes        0,MatTranspose_MPIDense,
80239ddd567SLois Curfman McInnes        MatGetInfo_MPIDense,0,
803096963f5SLois Curfman McInnes        MatGetDiagonal_MPIDense,0,MatNorm_MPIDense,
80439ddd567SLois Curfman McInnes        MatAssemblyBegin_MPIDense,MatAssemblyEnd_MPIDense,
8058965ea79SLois Curfman McInnes        0,
80639ddd567SLois Curfman McInnes        MatSetOption_MPIDense,MatZeroEntries_MPIDense,MatZeroRows_MPIDense,
8078c469469SLois Curfman McInnes        0,0,0,
8088c469469SLois Curfman McInnes /*       0,MatLUFactorSymbolic_MPIDense,MatLUFactorNumeric_MPIDense, */
8098aaee692SLois Curfman McInnes        0,0,
81039ddd567SLois Curfman McInnes        MatGetSize_MPIDense,MatGetLocalSize_MPIDense,
81139ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
812ff14e315SSatish Balay        0,0, MatGetArray_MPIDense, MatRestoreArray_MPIDense,
813905e6a2fSBarry Smith        0,MatConvertSameType_MPIDense,
814b49de8d1SLois Curfman McInnes        0,0,0,0,0,
815*4e220ebcSLois Curfman McInnes        0,0,MatGetValues_MPIDense,0,0,MatScale_MPIDense,
816*4e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_MPIDense};
8178965ea79SLois Curfman McInnes 
8188965ea79SLois Curfman McInnes /*@C
81939ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
8208965ea79SLois Curfman McInnes 
8218965ea79SLois Curfman McInnes    Input Parameters:
8228965ea79SLois Curfman McInnes .  comm - MPI communicator
8238965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
8248965ea79SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated
8258965ea79SLois Curfman McInnes            if N is given)
8268965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
8278965ea79SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated
8288965ea79SLois Curfman McInnes            if n is given)
829b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
830dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
8318965ea79SLois Curfman McInnes 
8328965ea79SLois Curfman McInnes    Output Parameter:
833477f1c0bSLois Curfman McInnes .  A - the matrix
8348965ea79SLois Curfman McInnes 
8358965ea79SLois Curfman McInnes    Notes:
83639ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
83739ddd567SLois Curfman McInnes    storage by columns.
8388965ea79SLois Curfman McInnes 
83918f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
84018f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
841b4fd4287SBarry Smith    set data=PETSC_NULL.
84218f449edSLois Curfman McInnes 
8438965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
8448965ea79SLois Curfman McInnes    (possibly both).
8458965ea79SLois Curfman McInnes 
8463501a2bdSLois Curfman McInnes    Currently, the only parallel dense matrix decomposition is by rows,
8473501a2bdSLois Curfman McInnes    so that n=N and each submatrix owns all of the global columns.
8483501a2bdSLois Curfman McInnes 
84939ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel
8508965ea79SLois Curfman McInnes 
85139ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
8528965ea79SLois Curfman McInnes @*/
853477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
8548965ea79SLois Curfman McInnes {
8558965ea79SLois Curfman McInnes   Mat          mat;
85639ddd567SLois Curfman McInnes   Mat_MPIDense *a;
85725cdf11fSBarry Smith   int          ierr, i,flg;
8588965ea79SLois Curfman McInnes 
859ed2daf61SLois Curfman McInnes   /* Note:  For now, when data is specified above, this assumes the user correctly
860ed2daf61SLois Curfman McInnes    allocates the local dense storage space.  We should add error checking. */
86118f449edSLois Curfman McInnes 
862477f1c0bSLois Curfman McInnes   *A = 0;
8630452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIDENSE,comm);
8648965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
8650452661fSBarry Smith   mat->data       = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
8668965ea79SLois Curfman McInnes   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
86739ddd567SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIDense;
86839ddd567SLois Curfman McInnes   mat->view       = MatView_MPIDense;
8698965ea79SLois Curfman McInnes   mat->factor     = 0;
8708965ea79SLois Curfman McInnes 
871622d7880SLois Curfman McInnes   a->factor     = 0;
8728965ea79SLois Curfman McInnes   a->insertmode = NOT_SET_VALUES;
8738965ea79SLois Curfman McInnes   MPI_Comm_rank(comm,&a->rank);
8748965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&a->size);
8758965ea79SLois Curfman McInnes 
87639ddd567SLois Curfman McInnes   if (M == PETSC_DECIDE) MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm);
8778965ea79SLois Curfman McInnes   if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);}
87839ddd567SLois Curfman McInnes 
87939ddd567SLois Curfman McInnes   /* each row stores all columns */
88039ddd567SLois Curfman McInnes   if (N == PETSC_DECIDE) N = n;
881d7e8b826SBarry Smith   if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);}
882d7e8b826SBarry Smith   /*  if (n != N) SETERRQ(1,"MatCreateMPIDense:For now, only n=N is supported"); */
883aca0ad90SLois Curfman McInnes   a->N = mat->N = N;
884aca0ad90SLois Curfman McInnes   a->M = mat->M = M;
885aca0ad90SLois Curfman McInnes   a->m = mat->m = m;
886aca0ad90SLois Curfman McInnes   a->n = mat->n = n;
8878965ea79SLois Curfman McInnes 
8888965ea79SLois Curfman McInnes   /* build local table of row and column ownerships */
889d7e8b826SBarry Smith   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners);
890d7e8b826SBarry Smith   a->cowners = a->rowners + a->size + 1;
8917056b6fcSBarry Smith   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIDense));
8928965ea79SLois Curfman McInnes   MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);
8938965ea79SLois Curfman McInnes   a->rowners[0] = 0;
8948965ea79SLois Curfman McInnes   for ( i=2; i<=a->size; i++ ) {
8958965ea79SLois Curfman McInnes     a->rowners[i] += a->rowners[i-1];
8968965ea79SLois Curfman McInnes   }
8978965ea79SLois Curfman McInnes   a->rstart = a->rowners[a->rank];
8988965ea79SLois Curfman McInnes   a->rend   = a->rowners[a->rank+1];
899d7e8b826SBarry Smith   MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);
900d7e8b826SBarry Smith   a->cowners[0] = 0;
901d7e8b826SBarry Smith   for ( i=2; i<=a->size; i++ ) {
902d7e8b826SBarry Smith     a->cowners[i] += a->cowners[i-1];
903d7e8b826SBarry Smith   }
9048965ea79SLois Curfman McInnes 
90518f449edSLois Curfman McInnes   ierr = MatCreateSeqDense(MPI_COMM_SELF,m,N,data,&a->A); CHKERRQ(ierr);
9068965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
9078965ea79SLois Curfman McInnes 
9088965ea79SLois Curfman McInnes   /* build cache for off array entries formed */
9098965ea79SLois Curfman McInnes   ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr);
9108965ea79SLois Curfman McInnes 
9118965ea79SLois Curfman McInnes   /* stuff used for matrix vector multiply */
9128965ea79SLois Curfman McInnes   a->lvec        = 0;
9138965ea79SLois Curfman McInnes   a->Mvctx       = 0;
91439b7565bSBarry Smith   a->roworiented = 1;
9158965ea79SLois Curfman McInnes 
916477f1c0bSLois Curfman McInnes   *A = mat;
91725cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
91825cdf11fSBarry Smith   if (flg) {
9198c469469SLois Curfman McInnes     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
9208c469469SLois Curfman McInnes   }
9218965ea79SLois Curfman McInnes   return 0;
9228965ea79SLois Curfman McInnes }
9238965ea79SLois Curfman McInnes 
9243d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat A,Mat *newmat,int cpvalues)
9258965ea79SLois Curfman McInnes {
9268965ea79SLois Curfman McInnes   Mat          mat;
9273501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data;
92839ddd567SLois Curfman McInnes   int          ierr;
9292ba99913SLois Curfman McInnes   FactorCtx    *factor;
9308965ea79SLois Curfman McInnes 
9318965ea79SLois Curfman McInnes   *newmat       = 0;
9320452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIDENSE,A->comm);
9338965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
9340452661fSBarry Smith   mat->data      = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
9358965ea79SLois Curfman McInnes   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
93639ddd567SLois Curfman McInnes   mat->destroy   = MatDestroy_MPIDense;
93739ddd567SLois Curfman McInnes   mat->view      = MatView_MPIDense;
9383501a2bdSLois Curfman McInnes   mat->factor    = A->factor;
939c456f294SBarry Smith   mat->assembled = PETSC_TRUE;
9408965ea79SLois Curfman McInnes 
94144cd7ae7SLois Curfman McInnes   a->m = mat->m = oldmat->m;
94244cd7ae7SLois Curfman McInnes   a->n = mat->n = oldmat->n;
94344cd7ae7SLois Curfman McInnes   a->M = mat->M = oldmat->M;
94444cd7ae7SLois Curfman McInnes   a->N = mat->N = oldmat->N;
9452ba99913SLois Curfman McInnes   if (oldmat->factor) {
9462ba99913SLois Curfman McInnes     a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx)); CHKPTRQ(factor);
9472ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
9482ba99913SLois Curfman McInnes   } else a->factor = 0;
9498965ea79SLois Curfman McInnes 
9508965ea79SLois Curfman McInnes   a->rstart     = oldmat->rstart;
9518965ea79SLois Curfman McInnes   a->rend       = oldmat->rend;
9528965ea79SLois Curfman McInnes   a->size       = oldmat->size;
9538965ea79SLois Curfman McInnes   a->rank       = oldmat->rank;
9548965ea79SLois Curfman McInnes   a->insertmode = NOT_SET_VALUES;
9558965ea79SLois Curfman McInnes 
9560452661fSBarry Smith   a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners);
95739ddd567SLois Curfman McInnes   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIDense));
9588965ea79SLois Curfman McInnes   PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));
9598965ea79SLois Curfman McInnes   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
9608965ea79SLois Curfman McInnes 
9618965ea79SLois Curfman McInnes   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
9628965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->lvec);
96355659b69SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
9648965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->Mvctx);
9658965ea79SLois Curfman McInnes   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
9668965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
9678965ea79SLois Curfman McInnes   *newmat = mat;
9688965ea79SLois Curfman McInnes   return 0;
9698965ea79SLois Curfman McInnes }
9708965ea79SLois Curfman McInnes 
97177c4ece6SBarry Smith #include "sys.h"
9728965ea79SLois Curfman McInnes 
97390ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat)
97490ace30eSBarry Smith {
97590ace30eSBarry Smith   int        *rowners, i,size,rank,m,rstart,rend,ierr,nz,j;
97690ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
97790ace30eSBarry Smith   MPI_Status status;
97890ace30eSBarry Smith 
97990ace30eSBarry Smith   MPI_Comm_rank(comm,&rank);
98090ace30eSBarry Smith   MPI_Comm_size(comm,&size);
98190ace30eSBarry Smith 
98290ace30eSBarry Smith   /* determine ownership of all rows */
98390ace30eSBarry Smith   m = M/size + ((M % size) > rank);
98490ace30eSBarry Smith   rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
98590ace30eSBarry Smith   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
98690ace30eSBarry Smith   rowners[0] = 0;
98790ace30eSBarry Smith   for ( i=2; i<=size; i++ ) {
98890ace30eSBarry Smith     rowners[i] += rowners[i-1];
98990ace30eSBarry Smith   }
99090ace30eSBarry Smith   rstart = rowners[rank];
99190ace30eSBarry Smith   rend   = rowners[rank+1];
99290ace30eSBarry Smith 
99390ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
99490ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array); CHKERRQ(ierr);
99590ace30eSBarry Smith 
99690ace30eSBarry Smith   if (!rank) {
99790ace30eSBarry Smith     vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
99890ace30eSBarry Smith 
99990ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
100077c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,BINARY_SCALAR); CHKERRQ(ierr);
100190ace30eSBarry Smith 
100290ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
100390ace30eSBarry Smith     vals_ptr = vals;
100490ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
100590ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
100690ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
100790ace30eSBarry Smith       }
100890ace30eSBarry Smith     }
100990ace30eSBarry Smith 
101090ace30eSBarry Smith     /* read in other processors and ship out */
101190ace30eSBarry Smith     for ( i=1; i<size; i++ ) {
101290ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
101377c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
101490ace30eSBarry Smith       MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);
101590ace30eSBarry Smith     }
101690ace30eSBarry Smith   }
101790ace30eSBarry Smith   else {
101890ace30eSBarry Smith     /* receive numeric values */
101990ace30eSBarry Smith     vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
102090ace30eSBarry Smith 
102190ace30eSBarry Smith     /* receive message of values*/
102290ace30eSBarry Smith     MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);
102390ace30eSBarry Smith 
102490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
102590ace30eSBarry Smith     vals_ptr = vals;
102690ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
102790ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
102890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
102990ace30eSBarry Smith       }
103090ace30eSBarry Smith     }
103190ace30eSBarry Smith   }
103290ace30eSBarry Smith   PetscFree(rowners);
103390ace30eSBarry Smith   PetscFree(vals);
10346d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10356d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
103690ace30eSBarry Smith   return 0;
103790ace30eSBarry Smith }
103890ace30eSBarry Smith 
103990ace30eSBarry Smith 
104019bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
10418965ea79SLois Curfman McInnes {
10428965ea79SLois Curfman McInnes   Mat          A;
10438965ea79SLois Curfman McInnes   int          i, nz, ierr, j,rstart, rend, fd;
10448965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
104519bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
10468965ea79SLois Curfman McInnes   MPI_Status   status;
10478965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
10488965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
104919bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
10508965ea79SLois Curfman McInnes 
10518965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank);
10528965ea79SLois Curfman McInnes   if (!rank) {
105390ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
105477c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr);
105539ddd567SLois Curfman McInnes     if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:not matrix object");
10568965ea79SLois Curfman McInnes   }
10578965ea79SLois Curfman McInnes 
10588965ea79SLois Curfman McInnes   MPI_Bcast(header+1,3,MPI_INT,0,comm);
105990ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
106090ace30eSBarry Smith 
106190ace30eSBarry Smith   /*
106290ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
106390ace30eSBarry Smith   */
106490ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
106590ace30eSBarry Smith     return MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);
106690ace30eSBarry Smith   }
106790ace30eSBarry Smith 
10688965ea79SLois Curfman McInnes   /* determine ownership of all rows */
10698965ea79SLois Curfman McInnes   m = M/size + ((M % size) > rank);
10700452661fSBarry Smith   rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
10718965ea79SLois Curfman McInnes   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
10728965ea79SLois Curfman McInnes   rowners[0] = 0;
10738965ea79SLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
10748965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
10758965ea79SLois Curfman McInnes   }
10768965ea79SLois Curfman McInnes   rstart = rowners[rank];
10778965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
10788965ea79SLois Curfman McInnes 
10798965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
10800452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
10818965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
10828965ea79SLois Curfman McInnes   if (!rank) {
10830452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
108477c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
10850452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts);
10868965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
10878965ea79SLois Curfman McInnes     MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
10880452661fSBarry Smith     PetscFree(sndcounts);
10898965ea79SLois Curfman McInnes   }
10908965ea79SLois Curfman McInnes   else {
10918965ea79SLois Curfman McInnes     MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);
10928965ea79SLois Curfman McInnes   }
10938965ea79SLois Curfman McInnes 
10948965ea79SLois Curfman McInnes   if (!rank) {
10958965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
10960452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz);
1097cddf8d76SBarry Smith     PetscMemzero(procsnz,size*sizeof(int));
10988965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
10998965ea79SLois Curfman McInnes       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
11008965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
11018965ea79SLois Curfman McInnes       }
11028965ea79SLois Curfman McInnes     }
11030452661fSBarry Smith     PetscFree(rowlengths);
11048965ea79SLois Curfman McInnes 
11058965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
11068965ea79SLois Curfman McInnes     maxnz = 0;
11078965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
11080452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
11098965ea79SLois Curfman McInnes     }
11100452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols);
11118965ea79SLois Curfman McInnes 
11128965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
11138965ea79SLois Curfman McInnes     nz = procsnz[0];
11140452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
111577c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr);
11168965ea79SLois Curfman McInnes 
11178965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
11188965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
11198965ea79SLois Curfman McInnes       nz = procsnz[i];
112077c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
11218965ea79SLois Curfman McInnes       MPI_Send(cols,nz,MPI_INT,i,tag,comm);
11228965ea79SLois Curfman McInnes     }
11230452661fSBarry Smith     PetscFree(cols);
11248965ea79SLois Curfman McInnes   }
11258965ea79SLois Curfman McInnes   else {
11268965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
11278965ea79SLois Curfman McInnes     nz = 0;
11288965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
11298965ea79SLois Curfman McInnes       nz += ourlens[i];
11308965ea79SLois Curfman McInnes     }
11310452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
11328965ea79SLois Curfman McInnes 
11338965ea79SLois Curfman McInnes     /* receive message of column indices*/
11348965ea79SLois Curfman McInnes     MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);
11358965ea79SLois Curfman McInnes     MPI_Get_count(&status,MPI_INT,&maxnz);
113639ddd567SLois Curfman McInnes     if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file");
11378965ea79SLois Curfman McInnes   }
11388965ea79SLois Curfman McInnes 
11398965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1140cddf8d76SBarry Smith   PetscMemzero(offlens,m*sizeof(int));
11418965ea79SLois Curfman McInnes   jj = 0;
11428965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
11438965ea79SLois Curfman McInnes     for ( j=0; j<ourlens[i]; j++ ) {
11448965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
11458965ea79SLois Curfman McInnes       jj++;
11468965ea79SLois Curfman McInnes     }
11478965ea79SLois Curfman McInnes   }
11488965ea79SLois Curfman McInnes 
11498965ea79SLois Curfman McInnes   /* create our matrix */
11508965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
11518965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
11528965ea79SLois Curfman McInnes   }
1153b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
11548965ea79SLois Curfman McInnes   A = *newmat;
11558965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
11568965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
11578965ea79SLois Curfman McInnes   }
11588965ea79SLois Curfman McInnes 
11598965ea79SLois Curfman McInnes   if (!rank) {
11600452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
11618965ea79SLois Curfman McInnes 
11628965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
11638965ea79SLois Curfman McInnes     nz = procsnz[0];
116477c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
11658965ea79SLois Curfman McInnes 
11668965ea79SLois Curfman McInnes     /* insert into matrix */
11678965ea79SLois Curfman McInnes     jj      = rstart;
11688965ea79SLois Curfman McInnes     smycols = mycols;
11698965ea79SLois Curfman McInnes     svals   = vals;
11708965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
11718965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
11728965ea79SLois Curfman McInnes       smycols += ourlens[i];
11738965ea79SLois Curfman McInnes       svals   += ourlens[i];
11748965ea79SLois Curfman McInnes       jj++;
11758965ea79SLois Curfman McInnes     }
11768965ea79SLois Curfman McInnes 
11778965ea79SLois Curfman McInnes     /* read in other processors and ship out */
11788965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
11798965ea79SLois Curfman McInnes       nz = procsnz[i];
118077c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
11818965ea79SLois Curfman McInnes       MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
11828965ea79SLois Curfman McInnes     }
11830452661fSBarry Smith     PetscFree(procsnz);
11848965ea79SLois Curfman McInnes   }
11858965ea79SLois Curfman McInnes   else {
11868965ea79SLois Curfman McInnes     /* receive numeric values */
11870452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
11888965ea79SLois Curfman McInnes 
11898965ea79SLois Curfman McInnes     /* receive message of values*/
11908965ea79SLois Curfman McInnes     MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
11918965ea79SLois Curfman McInnes     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
119239ddd567SLois Curfman McInnes     if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file");
11938965ea79SLois Curfman McInnes 
11948965ea79SLois Curfman McInnes     /* insert into matrix */
11958965ea79SLois Curfman McInnes     jj      = rstart;
11968965ea79SLois Curfman McInnes     smycols = mycols;
11978965ea79SLois Curfman McInnes     svals   = vals;
11988965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
11998965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
12008965ea79SLois Curfman McInnes       smycols += ourlens[i];
12018965ea79SLois Curfman McInnes       svals   += ourlens[i];
12028965ea79SLois Curfman McInnes       jj++;
12038965ea79SLois Curfman McInnes     }
12048965ea79SLois Curfman McInnes   }
12050452661fSBarry Smith   PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners);
12068965ea79SLois Curfman McInnes 
12076d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12086d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12098965ea79SLois Curfman McInnes   return 0;
12108965ea79SLois Curfman McInnes }
121190ace30eSBarry Smith 
121290ace30eSBarry Smith 
121390ace30eSBarry Smith 
121490ace30eSBarry Smith 
121590ace30eSBarry Smith 
1216