xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision e20fef11c1a1457fc77d0de77d14a2dd6b3370e4)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*e20fef11SSatish Balay static char vcid[] = "$Id: mpidense.c,v 1.91 1998/05/29 20:37:01 bsmith Exp balay $";
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 
125615d1e5SSatish Balay #undef __FUNC__
135615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIDense"
148f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
158965ea79SLois Curfman McInnes {
1639b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense *) mat->data;
1739b7565bSBarry Smith   int          ierr, i, j, rstart = A->rstart, rend = A->rend, row;
1839b7565bSBarry Smith   int          roworiented = A->roworiented;
198965ea79SLois Curfman McInnes 
203a40ed3dSBarry Smith   PetscFunctionBegin;
218965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
22a8c6a408SBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
23a8c6a408SBarry Smith     if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
248965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
258965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
2639b7565bSBarry Smith       if (roworiented) {
2739b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv); CHKERRQ(ierr);
283a40ed3dSBarry Smith       } else {
298965ea79SLois Curfman McInnes         for ( j=0; j<n; j++ ) {
30a8c6a408SBarry Smith           if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
31a8c6a408SBarry Smith           if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
3239b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv); CHKERRQ(ierr);
3339b7565bSBarry Smith         }
348965ea79SLois Curfman McInnes       }
353a40ed3dSBarry Smith     } else {
3639b7565bSBarry Smith       if (roworiented) {
3739b7565bSBarry Smith         ierr = StashValues_Private(&A->stash,idxm[i],n,idxn,v+i*n,addv); CHKERRQ(ierr);
383a40ed3dSBarry Smith       } else { /* must stash each seperately */
3939b7565bSBarry Smith         row = idxm[i];
4039b7565bSBarry Smith         for ( j=0; j<n; j++ ) {
417056b6fcSBarry Smith           ierr = StashValues_Private(&A->stash,row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
4239b7565bSBarry Smith         }
4339b7565bSBarry Smith       }
44b49de8d1SLois Curfman McInnes     }
45b49de8d1SLois Curfman McInnes   }
463a40ed3dSBarry Smith   PetscFunctionReturn(0);
47b49de8d1SLois Curfman McInnes }
48b49de8d1SLois Curfman McInnes 
495615d1e5SSatish Balay #undef __FUNC__
505615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIDense"
518f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
52b49de8d1SLois Curfman McInnes {
53b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
54b49de8d1SLois Curfman McInnes   int          ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row;
55b49de8d1SLois Curfman McInnes 
563a40ed3dSBarry Smith   PetscFunctionBegin;
57b49de8d1SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
58a8c6a408SBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
59a8c6a408SBarry Smith     if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
60b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
61b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
62b49de8d1SLois Curfman McInnes       for ( j=0; j<n; j++ ) {
63a8c6a408SBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
64a8c6a408SBarry Smith         if (idxn[j] >= mdn->N) {
65a8c6a408SBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
66a8c6a408SBarry Smith         }
67b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j); CHKERRQ(ierr);
68b49de8d1SLois Curfman McInnes       }
69a8c6a408SBarry Smith     } else {
70a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported");
718965ea79SLois Curfman McInnes     }
728965ea79SLois Curfman McInnes   }
733a40ed3dSBarry Smith   PetscFunctionReturn(0);
748965ea79SLois Curfman McInnes }
758965ea79SLois Curfman McInnes 
765615d1e5SSatish Balay #undef __FUNC__
775615d1e5SSatish Balay #define __FUNC__ "MatGetArray_MPIDense"
788f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array)
79ff14e315SSatish Balay {
80ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
81ff14e315SSatish Balay   int          ierr;
82ff14e315SSatish Balay 
833a40ed3dSBarry Smith   PetscFunctionBegin;
84ff14e315SSatish Balay   ierr = MatGetArray(a->A,array); CHKERRQ(ierr);
853a40ed3dSBarry Smith   PetscFunctionReturn(0);
86ff14e315SSatish Balay }
87ff14e315SSatish Balay 
885615d1e5SSatish Balay #undef __FUNC__
895615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_MPIDense"
908f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array)
91ff14e315SSatish Balay {
923a40ed3dSBarry Smith   PetscFunctionBegin;
933a40ed3dSBarry Smith   PetscFunctionReturn(0);
94ff14e315SSatish Balay }
95ff14e315SSatish Balay 
965615d1e5SSatish Balay #undef __FUNC__
975615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIDense"
988f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
998965ea79SLois Curfman McInnes {
10039ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
1018965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
10239ddd567SLois Curfman McInnes   int          size = mdn->size, *owners = mdn->rowners, rank = mdn->rank;
1038965ea79SLois Curfman McInnes   int          *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
10439ddd567SLois Curfman McInnes   int          tag = mat->tag, *owner,*starts,count,ierr;
1058965ea79SLois Curfman McInnes   InsertMode   addv;
10639ddd567SLois Curfman McInnes   MPI_Request  *send_waits,*recv_waits;
1078965ea79SLois Curfman McInnes   Scalar       *rvalues,*svalues;
1088965ea79SLois Curfman McInnes 
1093a40ed3dSBarry Smith   PetscFunctionBegin;
1108965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
111ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1127056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
113a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs");
1148965ea79SLois Curfman McInnes   }
115e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
1168965ea79SLois Curfman McInnes 
1178965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
1180452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
119cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
1200452661fSBarry Smith   owner = (int *) PetscMalloc( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
12139ddd567SLois Curfman McInnes   for ( i=0; i<mdn->stash.n; i++ ) {
12239ddd567SLois Curfman McInnes     idx = mdn->stash.idx[i];
1238965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
1248965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
1258965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
1268965ea79SLois Curfman McInnes       }
1278965ea79SLois Curfman McInnes     }
1288965ea79SLois Curfman McInnes   }
1298965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
1308965ea79SLois Curfman McInnes 
1318965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
1320452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
133ca161407SBarry Smith   ierr = MPI_Allreduce(procs,work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
1348965ea79SLois Curfman McInnes   nreceives = work[rank];
135a8c6a408SBarry Smith   if (nreceives > size) SETERRQ(PETSC_ERR_PLIB,0,"Internal PETSc error");
136ca161407SBarry Smith   ierr = MPI_Allreduce(nprocs,work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
1378965ea79SLois Curfman McInnes   nmax = work[rank];
1380452661fSBarry Smith   PetscFree(work);
1398965ea79SLois Curfman McInnes 
1408965ea79SLois Curfman McInnes   /* post receives:
1418965ea79SLois Curfman McInnes        1) each message will consist of ordered pairs
1428965ea79SLois Curfman McInnes      (global index,value) we store the global index as a double
1438965ea79SLois Curfman McInnes      to simplify the message passing.
1448965ea79SLois Curfman McInnes        2) since we don't know how long each individual message is we
1458965ea79SLois Curfman McInnes      allocate the largest needed buffer for each receive. Potentially
1468965ea79SLois Curfman McInnes      this is a lot of wasted space.
1478965ea79SLois Curfman McInnes 
1488965ea79SLois Curfman McInnes        This could be done better.
1498965ea79SLois Curfman McInnes   */
1503b2fbd54SBarry Smith   rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));CHKPTRQ(rvalues);
1513b2fbd54SBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
1528965ea79SLois Curfman McInnes   for ( i=0; i<nreceives; i++ ) {
153ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
1548965ea79SLois Curfman McInnes   }
1558965ea79SLois Curfman McInnes 
1568965ea79SLois Curfman McInnes   /* do sends:
1578965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
1588965ea79SLois Curfman McInnes          the ith processor
1598965ea79SLois Curfman McInnes   */
1603b2fbd54SBarry Smith   svalues = (Scalar *) PetscMalloc( 3*(mdn->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues);
1613b2fbd54SBarry Smith   send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
1620452661fSBarry Smith   starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts);
1638965ea79SLois Curfman McInnes   starts[0] = 0;
1648965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
16539ddd567SLois Curfman McInnes   for ( i=0; i<mdn->stash.n; i++ ) {
16639ddd567SLois Curfman McInnes     svalues[3*starts[owner[i]]]       = (Scalar)  mdn->stash.idx[i];
16739ddd567SLois Curfman McInnes     svalues[3*starts[owner[i]]+1]     = (Scalar)  mdn->stash.idy[i];
16839ddd567SLois Curfman McInnes     svalues[3*(starts[owner[i]]++)+2] =  mdn->stash.array[i];
1698965ea79SLois Curfman McInnes   }
1700452661fSBarry Smith   PetscFree(owner);
1718965ea79SLois Curfman McInnes   starts[0] = 0;
1728965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1738965ea79SLois Curfman McInnes   count = 0;
1748965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
1758965ea79SLois Curfman McInnes     if (procs[i]) {
176ca161407SBarry Smith       ierr = MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
1778965ea79SLois Curfman McInnes     }
1788965ea79SLois Curfman McInnes   }
1790452661fSBarry Smith   PetscFree(starts); PetscFree(nprocs);
1808965ea79SLois Curfman McInnes 
1818965ea79SLois Curfman McInnes   /* Free cache space */
182d2dc9b81SLois Curfman McInnes   PLogInfo(mat,"MatAssemblyBegin_MPIDense:Number of off-processor values %d\n",mdn->stash.n);
18339ddd567SLois Curfman McInnes   ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr);
1848965ea79SLois Curfman McInnes 
18539ddd567SLois Curfman McInnes   mdn->svalues    = svalues;    mdn->rvalues = rvalues;
18639ddd567SLois Curfman McInnes   mdn->nsends     = nsends;     mdn->nrecvs = nreceives;
18739ddd567SLois Curfman McInnes   mdn->send_waits = send_waits; mdn->recv_waits = recv_waits;
18839ddd567SLois Curfman McInnes   mdn->rmax       = nmax;
1898965ea79SLois Curfman McInnes 
1903a40ed3dSBarry Smith   PetscFunctionReturn(0);
1918965ea79SLois Curfman McInnes }
19239ddd567SLois Curfman McInnes extern int MatSetUpMultiply_MPIDense(Mat);
1938965ea79SLois Curfman McInnes 
1945615d1e5SSatish Balay #undef __FUNC__
1955615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIDense"
1968f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
1978965ea79SLois Curfman McInnes {
19839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
1998965ea79SLois Curfman McInnes   MPI_Status   *send_status,recv_status;
20039ddd567SLois Curfman McInnes   int          imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col;
2018965ea79SLois Curfman McInnes   Scalar       *values,val;
202e0fa3b82SLois Curfman McInnes   InsertMode   addv = mat->insertmode;
2038965ea79SLois Curfman McInnes 
2043a40ed3dSBarry Smith   PetscFunctionBegin;
2058965ea79SLois Curfman McInnes   /*  wait on receives */
2068965ea79SLois Curfman McInnes   while (count) {
207ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
2088965ea79SLois Curfman McInnes     /* unpack receives into our local space */
20939ddd567SLois Curfman McInnes     values = mdn->rvalues + 3*imdex*mdn->rmax;
210ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_SCALAR,&n);CHKERRQ(ierr);
2118965ea79SLois Curfman McInnes     n = n/3;
2128965ea79SLois Curfman McInnes     for ( i=0; i<n; i++ ) {
213227d817aSBarry Smith       row = (int) PetscReal(values[3*i]) - mdn->rstart;
214227d817aSBarry Smith       col = (int) PetscReal(values[3*i+1]);
2158965ea79SLois Curfman McInnes       val = values[3*i+2];
21639ddd567SLois Curfman McInnes       if (col >= 0 && col < mdn->N) {
21739ddd567SLois Curfman McInnes         MatSetValues(mdn->A,1,&row,1,&col,&val,addv);
2188965ea79SLois Curfman McInnes       }
219a8c6a408SBarry Smith       else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid column");}
2208965ea79SLois Curfman McInnes     }
2218965ea79SLois Curfman McInnes     count--;
2228965ea79SLois Curfman McInnes   }
2230452661fSBarry Smith   PetscFree(mdn->recv_waits); PetscFree(mdn->rvalues);
2248965ea79SLois Curfman McInnes 
2258965ea79SLois Curfman McInnes   /* wait on sends */
22639ddd567SLois Curfman McInnes   if (mdn->nsends) {
2277056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(mdn->nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
228ca161407SBarry Smith     ierr        = MPI_Waitall(mdn->nsends,mdn->send_waits,send_status);CHKERRQ(ierr);
2290452661fSBarry Smith     PetscFree(send_status);
2308965ea79SLois Curfman McInnes   }
2310452661fSBarry Smith   PetscFree(mdn->send_waits); PetscFree(mdn->svalues);
2328965ea79SLois Curfman McInnes 
23339ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr);
23439ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr);
2358965ea79SLois Curfman McInnes 
2366d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23739ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr);
2388965ea79SLois Curfman McInnes   }
2393a40ed3dSBarry Smith   PetscFunctionReturn(0);
2408965ea79SLois Curfman McInnes }
2418965ea79SLois Curfman McInnes 
2425615d1e5SSatish Balay #undef __FUNC__
2435615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIDense"
2448f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2458965ea79SLois Curfman McInnes {
2463a40ed3dSBarry Smith   int          ierr;
24739ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense *) A->data;
2483a40ed3dSBarry Smith 
2493a40ed3dSBarry Smith   PetscFunctionBegin;
2503a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
2528965ea79SLois Curfman McInnes }
2538965ea79SLois Curfman McInnes 
2545615d1e5SSatish Balay #undef __FUNC__
2555615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_MPIDense"
2568f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2574e220ebcSLois Curfman McInnes {
2583a40ed3dSBarry Smith   PetscFunctionBegin;
2594e220ebcSLois Curfman McInnes   *bs = 1;
2603a40ed3dSBarry Smith   PetscFunctionReturn(0);
2614e220ebcSLois Curfman McInnes }
2624e220ebcSLois Curfman McInnes 
2638965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2648965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2658965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2663501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2678965ea79SLois Curfman McInnes    routine.
2688965ea79SLois Curfman McInnes */
2695615d1e5SSatish Balay #undef __FUNC__
2705615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIDense"
2718f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag)
2728965ea79SLois Curfman McInnes {
27339ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense *) A->data;
2748965ea79SLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
2758965ea79SLois Curfman McInnes   int            *procs,*nprocs,j,found,idx,nsends,*work;
2768965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2778965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2788965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2798965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2808965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2818965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2828965ea79SLois Curfman McInnes   IS             istmp;
2838965ea79SLois Curfman McInnes 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
28577c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
2868965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
2878965ea79SLois Curfman McInnes 
2888965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
2890452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
290cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
2910452661fSBarry Smith   owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
2928965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
2938965ea79SLois Curfman McInnes     idx = rows[i];
2948965ea79SLois Curfman McInnes     found = 0;
2958965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
2968965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
2978965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
2988965ea79SLois Curfman McInnes       }
2998965ea79SLois Curfman McInnes     }
300a8c6a408SBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range");
3018965ea79SLois Curfman McInnes   }
3028965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
3038965ea79SLois Curfman McInnes 
3048965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
3050452661fSBarry Smith   work   = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
306ca161407SBarry Smith   ierr   = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
3078965ea79SLois Curfman McInnes   nrecvs = work[rank];
308ca161407SBarry Smith   ierr   = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
3098965ea79SLois Curfman McInnes   nmax   = work[rank];
3100452661fSBarry Smith   PetscFree(work);
3118965ea79SLois Curfman McInnes 
3128965ea79SLois Curfman McInnes   /* post receives:   */
3133a40ed3dSBarry Smith   rvalues    = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues);
3143a40ed3dSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
3158965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
316ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3178965ea79SLois Curfman McInnes   }
3188965ea79SLois Curfman McInnes 
3198965ea79SLois Curfman McInnes   /* do sends:
3208965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3218965ea79SLois Curfman McInnes          the ith processor
3228965ea79SLois Curfman McInnes   */
3230452661fSBarry Smith   svalues    = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
3247056b6fcSBarry Smith   send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
3250452661fSBarry Smith   starts     = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts);
3268965ea79SLois Curfman McInnes   starts[0]  = 0;
3278965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3288965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
3298965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3308965ea79SLois Curfman McInnes   }
3318965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3328965ea79SLois Curfman McInnes 
3338965ea79SLois Curfman McInnes   starts[0] = 0;
3348965ea79SLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3358965ea79SLois Curfman McInnes   count = 0;
3368965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3378965ea79SLois Curfman McInnes     if (procs[i]) {
338ca161407SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3398965ea79SLois Curfman McInnes     }
3408965ea79SLois Curfman McInnes   }
3410452661fSBarry Smith   PetscFree(starts);
3428965ea79SLois Curfman McInnes 
3438965ea79SLois Curfman McInnes   base = owners[rank];
3448965ea79SLois Curfman McInnes 
3458965ea79SLois Curfman McInnes   /*  wait on receives */
3460452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
3478965ea79SLois Curfman McInnes   source = lens + nrecvs;
3488965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3498965ea79SLois Curfman McInnes   while (count) {
350ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3518965ea79SLois Curfman McInnes     /* unpack receives into our local space */
352ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3538965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3548965ea79SLois Curfman McInnes     lens[imdex]  = n;
3558965ea79SLois Curfman McInnes     slen += n;
3568965ea79SLois Curfman McInnes     count--;
3578965ea79SLois Curfman McInnes   }
3580452661fSBarry Smith   PetscFree(recv_waits);
3598965ea79SLois Curfman McInnes 
3608965ea79SLois Curfman McInnes   /* move the data into the send scatter */
3610452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
3628965ea79SLois Curfman McInnes   count = 0;
3638965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
3648965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3658965ea79SLois Curfman McInnes     for ( j=0; j<lens[i]; j++ ) {
3668965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3678965ea79SLois Curfman McInnes     }
3688965ea79SLois Curfman McInnes   }
3690452661fSBarry Smith   PetscFree(rvalues); PetscFree(lens);
3700452661fSBarry Smith   PetscFree(owner); PetscFree(nprocs);
3718965ea79SLois Curfman McInnes 
3728965ea79SLois Curfman McInnes   /* actually zap the local rows */
373029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
3748965ea79SLois Curfman McInnes   PLogObjectParent(A,istmp);
3750452661fSBarry Smith   PetscFree(lrows);
3768965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp); CHKERRQ(ierr);
3788965ea79SLois Curfman McInnes 
3798965ea79SLois Curfman McInnes   /* wait on sends */
3808965ea79SLois Curfman McInnes   if (nsends) {
3817056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
382ca161407SBarry Smith     ierr        = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
3830452661fSBarry Smith     PetscFree(send_status);
3848965ea79SLois Curfman McInnes   }
3850452661fSBarry Smith   PetscFree(send_waits); PetscFree(svalues);
3868965ea79SLois Curfman McInnes 
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
3888965ea79SLois Curfman McInnes }
3898965ea79SLois Curfman McInnes 
3905615d1e5SSatish Balay #undef __FUNC__
3915615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense"
3928f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3938965ea79SLois Curfman McInnes {
39439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
3958965ea79SLois Curfman McInnes   int          ierr;
396c456f294SBarry Smith 
3973a40ed3dSBarry Smith   PetscFunctionBegin;
39843a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39943a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40044cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy); CHKERRQ(ierr);
4013a40ed3dSBarry Smith   PetscFunctionReturn(0);
4028965ea79SLois Curfman McInnes }
4038965ea79SLois Curfman McInnes 
4045615d1e5SSatish Balay #undef __FUNC__
4055615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense"
4068f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4078965ea79SLois Curfman McInnes {
40839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4098965ea79SLois Curfman McInnes   int          ierr;
410c456f294SBarry Smith 
4113a40ed3dSBarry Smith   PetscFunctionBegin;
41243a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41343a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
41444cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr);
4153a40ed3dSBarry Smith   PetscFunctionReturn(0);
4168965ea79SLois Curfman McInnes }
4178965ea79SLois Curfman McInnes 
4185615d1e5SSatish Balay #undef __FUNC__
4195615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense"
4208f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy)
421096963f5SLois Curfman McInnes {
422096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
423096963f5SLois Curfman McInnes   int          ierr;
4243501a2bdSLois Curfman McInnes   Scalar       zero = 0.0;
425096963f5SLois Curfman McInnes 
4263a40ed3dSBarry Smith   PetscFunctionBegin;
4273501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy); CHKERRQ(ierr);
42844cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr);
429537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
430537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
4313a40ed3dSBarry Smith   PetscFunctionReturn(0);
432096963f5SLois Curfman McInnes }
433096963f5SLois Curfman McInnes 
4345615d1e5SSatish Balay #undef __FUNC__
4355615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense"
4368f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
437096963f5SLois Curfman McInnes {
438096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
439096963f5SLois Curfman McInnes   int          ierr;
440096963f5SLois Curfman McInnes 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
4423501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz); CHKERRQ(ierr);
44344cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr);
444537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
445537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
4463a40ed3dSBarry Smith   PetscFunctionReturn(0);
447096963f5SLois Curfman McInnes }
448096963f5SLois Curfman McInnes 
4495615d1e5SSatish Balay #undef __FUNC__
4505615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense"
4518f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4528965ea79SLois Curfman McInnes {
45339ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
454096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data;
45544cd7ae7SLois Curfman McInnes   int          ierr, len, i, n, m = a->m, radd;
45644cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
457ed3cc1f0SBarry Smith 
4583a40ed3dSBarry Smith   PetscFunctionBegin;
45944cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
460096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x); CHKERRQ(ierr);
461096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n); CHKERRQ(ierr);
462a8c6a408SBarry Smith   if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec");
46344cd7ae7SLois Curfman McInnes   len = PetscMin(aloc->m,aloc->n);
4647ddc982cSLois Curfman McInnes   radd = a->rstart*m;
46544cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
466096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
467096963f5SLois Curfman McInnes   }
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
4698965ea79SLois Curfman McInnes }
4708965ea79SLois Curfman McInnes 
4715615d1e5SSatish Balay #undef __FUNC__
4725615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense"
473e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4748965ea79SLois Curfman McInnes {
4753501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4768965ea79SLois Curfman McInnes   int          ierr;
477ed3cc1f0SBarry Smith 
4783a40ed3dSBarry Smith   PetscFunctionBegin;
47994d884c6SBarry Smith   if (--mat->refct > 0) PetscFunctionReturn(0);
48094d884c6SBarry Smith 
48194d884c6SBarry Smith   if (mat->mapping) {
48294d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
48394d884c6SBarry Smith   }
48494d884c6SBarry Smith   if (mat->bmapping) {
48594d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr);
48694d884c6SBarry Smith   }
4873a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
488e1311b90SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N);
4898965ea79SLois Curfman McInnes #endif
4900452661fSBarry Smith   PetscFree(mdn->rowners);
4913501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A); CHKERRQ(ierr);
4923501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4933501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
494622d7880SLois Curfman McInnes   if (mdn->factor) {
495622d7880SLois Curfman McInnes     if (mdn->factor->temp)   PetscFree(mdn->factor->temp);
496622d7880SLois Curfman McInnes     if (mdn->factor->tag)    PetscFree(mdn->factor->tag);
497622d7880SLois Curfman McInnes     if (mdn->factor->pivots) PetscFree(mdn->factor->pivots);
498622d7880SLois Curfman McInnes     PetscFree(mdn->factor);
499622d7880SLois Curfman McInnes   }
5000452661fSBarry Smith   PetscFree(mdn);
50190f02eecSBarry Smith   if (mat->mapping) {
50290f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
50390f02eecSBarry Smith   }
5048965ea79SLois Curfman McInnes   PLogObjectDestroy(mat);
5050452661fSBarry Smith   PetscHeaderDestroy(mat);
5063a40ed3dSBarry Smith   PetscFunctionReturn(0);
5078965ea79SLois Curfman McInnes }
50839ddd567SLois Curfman McInnes 
5098965ea79SLois Curfman McInnes #include "pinclude/pviewer.h"
5108965ea79SLois Curfman McInnes 
5115615d1e5SSatish Balay #undef __FUNC__
5125615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary"
51339ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
5148965ea79SLois Curfman McInnes {
51539ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
5168965ea79SLois Curfman McInnes   int          ierr;
5177056b6fcSBarry Smith 
5183a40ed3dSBarry Smith   PetscFunctionBegin;
51939ddd567SLois Curfman McInnes   if (mdn->size == 1) {
52039ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5218965ea79SLois Curfman McInnes   }
522a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported");
5233a40ed3dSBarry Smith   PetscFunctionReturn(0);
5248965ea79SLois Curfman McInnes }
5258965ea79SLois Curfman McInnes 
5265615d1e5SSatish Balay #undef __FUNC__
5275615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_ASCII"
52839ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer)
5298965ea79SLois Curfman McInnes {
53039ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
53139ddd567SLois Curfman McInnes   int          ierr, format;
5328965ea79SLois Curfman McInnes   FILE         *fd;
53319bcc07fSBarry Smith   ViewerType   vtype;
5348965ea79SLois Curfman McInnes 
5353a40ed3dSBarry Smith   PetscFunctionBegin;
5363a40ed3dSBarry Smith   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
53790ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
53890ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
539639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5404e220ebcSLois Curfman McInnes     int rank;
5414e220ebcSLois Curfman McInnes     MatInfo info;
542096963f5SLois Curfman McInnes     MPI_Comm_rank(mat->comm,&rank);
5434e220ebcSLois Curfman McInnes     ierr = MatGetInfo(mat,MAT_LOCAL,&info);
54477c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
5454e220ebcSLois Curfman McInnes       fprintf(fd,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m,
5464e220ebcSLois Curfman McInnes          (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
547096963f5SLois Curfman McInnes       fflush(fd);
54877c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
5493501a2bdSLois Curfman McInnes     ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr);
5503a40ed3dSBarry Smith     PetscFunctionReturn(0);
5513501a2bdSLois Curfman McInnes   }
55202cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_INFO) {
5533a40ed3dSBarry Smith     PetscFunctionReturn(0);
5548965ea79SLois Curfman McInnes   }
55519bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
55677c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
55739ddd567SLois Curfman McInnes     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n",
55839ddd567SLois Curfman McInnes              mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n);
55939ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5608965ea79SLois Curfman McInnes     fflush(fd);
56177c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
5623a40ed3dSBarry Smith   } else {
56339ddd567SLois Curfman McInnes     int size = mdn->size, rank = mdn->rank;
5648965ea79SLois Curfman McInnes     if (size == 1) {
56539ddd567SLois Curfman McInnes       ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5663a40ed3dSBarry Smith     } else {
5678965ea79SLois Curfman McInnes       /* assemble the entire matrix onto first processor. */
5688965ea79SLois Curfman McInnes       Mat          A;
56939ddd567SLois Curfman McInnes       int          M = mdn->M, N = mdn->N,m,row,i, nz, *cols;
57039ddd567SLois Curfman McInnes       Scalar       *vals;
57139ddd567SLois Curfman McInnes       Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data;
5728965ea79SLois Curfman McInnes 
5738965ea79SLois Curfman McInnes       if (!rank) {
5740513a670SBarry Smith         ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr);
5753a40ed3dSBarry Smith       } else {
5760513a670SBarry Smith         ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr);
5778965ea79SLois Curfman McInnes       }
5788965ea79SLois Curfman McInnes       PLogObjectParent(mat,A);
5798965ea79SLois Curfman McInnes 
58039ddd567SLois Curfman McInnes       /* Copy the matrix ... This isn't the most efficient means,
58139ddd567SLois Curfman McInnes          but it's quick for now */
58239ddd567SLois Curfman McInnes       row = mdn->rstart; m = Amdn->m;
5838965ea79SLois Curfman McInnes       for ( i=0; i<m; i++ ) {
58439ddd567SLois Curfman McInnes         ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
58539ddd567SLois Curfman McInnes         ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr);
58639ddd567SLois Curfman McInnes         ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
58739ddd567SLois Curfman McInnes         row++;
5888965ea79SLois Curfman McInnes       }
5898965ea79SLois Curfman McInnes 
5906d4a8577SBarry Smith       ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5916d4a8577SBarry Smith       ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5928965ea79SLois Curfman McInnes       if (!rank) {
59339ddd567SLois Curfman McInnes         ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr);
5948965ea79SLois Curfman McInnes       }
5958965ea79SLois Curfman McInnes       ierr = MatDestroy(A); CHKERRQ(ierr);
5968965ea79SLois Curfman McInnes     }
5978965ea79SLois Curfman McInnes   }
5983a40ed3dSBarry Smith   PetscFunctionReturn(0);
5998965ea79SLois Curfman McInnes }
6008965ea79SLois Curfman McInnes 
6015615d1e5SSatish Balay #undef __FUNC__
6025615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense"
603e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer)
6048965ea79SLois Curfman McInnes {
60539ddd567SLois Curfman McInnes   int          ierr;
606bcd2baecSBarry Smith   ViewerType   vtype;
6078965ea79SLois Curfman McInnes 
608bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
609bcd2baecSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
61039ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
6113a40ed3dSBarry Smith   } else if (vtype == ASCII_FILES_VIEWER) {
61239ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
6133a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
6143a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6155cd90555SBarry Smith   } else {
6165cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
6178965ea79SLois Curfman McInnes   }
6183a40ed3dSBarry Smith   PetscFunctionReturn(0);
6198965ea79SLois Curfman McInnes }
6208965ea79SLois Curfman McInnes 
6215615d1e5SSatish Balay #undef __FUNC__
6225615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense"
6238f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6248965ea79SLois Curfman McInnes {
6253501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6263501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6274e220ebcSLois Curfman McInnes   int          ierr;
6284e220ebcSLois Curfman McInnes   double       isend[5], irecv[5];
6298965ea79SLois Curfman McInnes 
6303a40ed3dSBarry Smith   PetscFunctionBegin;
6314e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
6324e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
6334e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
6344e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
6354e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6364e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info); CHKERRQ(ierr);
6374e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6384e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6398965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6404e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6414e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6424e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6434e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6444e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6458965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
646f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
6474e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6484e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6494e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6504e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6514e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6528965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
653f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
6544e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6554e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6564e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6574e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6584e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6598965ea79SLois Curfman McInnes   }
6604e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6614e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6624e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
6648965ea79SLois Curfman McInnes }
6658965ea79SLois Curfman McInnes 
6668c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*);
6678aaee692SLois Curfman McInnes    extern int MatLUFactorNumeric_MPIDense(Mat,Mat*);
6688aaee692SLois Curfman McInnes    extern int MatLUFactor_MPIDense(Mat,IS,IS,double);
6698aaee692SLois Curfman McInnes    extern int MatSolve_MPIDense(Mat,Vec,Vec);
6708c469469SLois Curfman McInnes    extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec);
6718aaee692SLois Curfman McInnes    extern int MatSolveTrans_MPIDense(Mat,Vec,Vec);
6728aaee692SLois Curfman McInnes    extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */
6738aaee692SLois Curfman McInnes 
6745615d1e5SSatish Balay #undef __FUNC__
6755615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense"
6768f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6778965ea79SLois Curfman McInnes {
67839ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
6798965ea79SLois Curfman McInnes 
6803a40ed3dSBarry Smith   PetscFunctionBegin;
6816d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6826d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
683c2653b3dSLois Curfman McInnes       op == MAT_NEW_NONZERO_LOCATION_ERROR ||
68496854ed6SLois Curfman McInnes       op == MAT_NEW_NONZERO_ALLOCATION_ERROR ||
685219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
686219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED) {
687b1fbbac0SLois Curfman McInnes         MatSetOption(a->A,op);
688b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
689aeafbbfcSLois Curfman McInnes         a->roworiented = 1;
6908965ea79SLois Curfman McInnes         MatSetOption(a->A,op);
691b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
692219d9a1aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
6936d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
6946d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
695b51ba29fSSatish Balay              op == MAT_YES_NEW_DIAGONALS ||
696b51ba29fSSatish Balay              op == MAT_USE_HASH_TABLE) {
697981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
6983a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
6993a40ed3dSBarry Smith     a->roworiented = 0; MatSetOption(a->A,op);
7003a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
7013a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7023a40ed3dSBarry Smith   } else {
7033a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7043a40ed3dSBarry Smith   }
7053a40ed3dSBarry Smith   PetscFunctionReturn(0);
7068965ea79SLois Curfman McInnes }
7078965ea79SLois Curfman McInnes 
7085615d1e5SSatish Balay #undef __FUNC__
7095615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense"
7108f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n)
7118965ea79SLois Curfman McInnes {
7123501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7133a40ed3dSBarry Smith 
7143a40ed3dSBarry Smith   PetscFunctionBegin;
7158965ea79SLois Curfman McInnes   *m = mat->M; *n = mat->N;
7163a40ed3dSBarry Smith   PetscFunctionReturn(0);
7178965ea79SLois Curfman McInnes }
7188965ea79SLois Curfman McInnes 
7195615d1e5SSatish Balay #undef __FUNC__
7205615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense"
7218f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n)
7228965ea79SLois Curfman McInnes {
7233501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7243a40ed3dSBarry Smith 
7253a40ed3dSBarry Smith   PetscFunctionBegin;
7268965ea79SLois Curfman McInnes   *m = mat->m; *n = mat->N;
7273a40ed3dSBarry Smith   PetscFunctionReturn(0);
7288965ea79SLois Curfman McInnes }
7298965ea79SLois Curfman McInnes 
7305615d1e5SSatish Balay #undef __FUNC__
7315615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense"
7328f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
7338965ea79SLois Curfman McInnes {
7343501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7353a40ed3dSBarry Smith 
7363a40ed3dSBarry Smith   PetscFunctionBegin;
7378965ea79SLois Curfman McInnes   *m = mat->rstart; *n = mat->rend;
7383a40ed3dSBarry Smith   PetscFunctionReturn(0);
7398965ea79SLois Curfman McInnes }
7408965ea79SLois Curfman McInnes 
7415615d1e5SSatish Balay #undef __FUNC__
7425615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense"
7438f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
7448965ea79SLois Curfman McInnes {
7453501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7463a40ed3dSBarry Smith   int          lrow, rstart = mat->rstart, rend = mat->rend,ierr;
7478965ea79SLois Curfman McInnes 
7483a40ed3dSBarry Smith   PetscFunctionBegin;
749a8c6a408SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows")
7508965ea79SLois Curfman McInnes   lrow = row - rstart;
7513a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7523a40ed3dSBarry Smith   PetscFunctionReturn(0);
7538965ea79SLois Curfman McInnes }
7548965ea79SLois Curfman McInnes 
7555615d1e5SSatish Balay #undef __FUNC__
7565615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense"
7578f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
7588965ea79SLois Curfman McInnes {
7593a40ed3dSBarry Smith   PetscFunctionBegin;
7600452661fSBarry Smith   if (idx) PetscFree(*idx);
7610452661fSBarry Smith   if (v) PetscFree(*v);
7623a40ed3dSBarry Smith   PetscFunctionReturn(0);
7638965ea79SLois Curfman McInnes }
7648965ea79SLois Curfman McInnes 
7655615d1e5SSatish Balay #undef __FUNC__
7665615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense"
7678f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm)
768096963f5SLois Curfman McInnes {
7693501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
7703501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
7713501a2bdSLois Curfman McInnes   int          ierr, i, j;
7723501a2bdSLois Curfman McInnes   double       sum = 0.0;
7733501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
7743501a2bdSLois Curfman McInnes 
7753a40ed3dSBarry Smith   PetscFunctionBegin;
7763501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
7773501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm); CHKERRQ(ierr);
7783501a2bdSLois Curfman McInnes   } else {
7793501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
7803501a2bdSLois Curfman McInnes       for (i=0; i<mat->n*mat->m; i++ ) {
7813a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
782*e20fef11SSatish Balay         sum += PetscReal(PetscConj(*v)*(*v)); v++;
7833501a2bdSLois Curfman McInnes #else
7843501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7853501a2bdSLois Curfman McInnes #endif
7863501a2bdSLois Curfman McInnes       }
787ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
7883501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
7893501a2bdSLois Curfman McInnes       PLogFlops(2*mat->n*mat->m);
7903a40ed3dSBarry Smith     } else if (type == NORM_1) {
7913501a2bdSLois Curfman McInnes       double *tmp, *tmp2;
7920452661fSBarry Smith       tmp  = (double *) PetscMalloc( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp);
7933501a2bdSLois Curfman McInnes       tmp2 = tmp + mdn->N;
794cddf8d76SBarry Smith       PetscMemzero(tmp,2*mdn->N*sizeof(double));
795096963f5SLois Curfman McInnes       *norm = 0.0;
7963501a2bdSLois Curfman McInnes       v = mat->v;
7973501a2bdSLois Curfman McInnes       for ( j=0; j<mat->n; j++ ) {
7983501a2bdSLois Curfman McInnes         for ( i=0; i<mat->m; i++ ) {
79967e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8003501a2bdSLois Curfman McInnes         }
8013501a2bdSLois Curfman McInnes       }
802ca161407SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
8033501a2bdSLois Curfman McInnes       for ( j=0; j<mdn->N; j++ ) {
8043501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
8053501a2bdSLois Curfman McInnes       }
8060452661fSBarry Smith       PetscFree(tmp);
8073501a2bdSLois Curfman McInnes       PLogFlops(mat->n*mat->m);
8083a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
8093501a2bdSLois Curfman McInnes       double ntemp;
8103501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr);
811ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
8123a40ed3dSBarry Smith     } else {
813a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
8143501a2bdSLois Curfman McInnes     }
8153501a2bdSLois Curfman McInnes   }
8163a40ed3dSBarry Smith   PetscFunctionReturn(0);
8173501a2bdSLois Curfman McInnes }
8183501a2bdSLois Curfman McInnes 
8195615d1e5SSatish Balay #undef __FUNC__
8205615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense"
8218f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8223501a2bdSLois Curfman McInnes {
8233501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
8243501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data;
8253501a2bdSLois Curfman McInnes   Mat          B;
8263501a2bdSLois Curfman McInnes   int          M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart;
8273501a2bdSLois Curfman McInnes   int          j, i, ierr;
8283501a2bdSLois Curfman McInnes   Scalar       *v;
8293501a2bdSLois Curfman McInnes 
8303a40ed3dSBarry Smith   PetscFunctionBegin;
8317056b6fcSBarry Smith   if (matout == PETSC_NULL && M != N) {
832a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place");
8337056b6fcSBarry Smith   }
8347056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8353501a2bdSLois Curfman McInnes 
8363501a2bdSLois Curfman McInnes   m = Aloc->m; n = Aloc->n; v = Aloc->v;
8370452661fSBarry Smith   rwork = (int *) PetscMalloc(n*sizeof(int)); CHKPTRQ(rwork);
8383501a2bdSLois Curfman McInnes   for ( j=0; j<n; j++ ) {
8393501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8403501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr);
8413501a2bdSLois Curfman McInnes     v   += m;
8423501a2bdSLois Curfman McInnes   }
8430452661fSBarry Smith   PetscFree(rwork);
8446d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8456d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8463638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
8473501a2bdSLois Curfman McInnes     *matout = B;
8483501a2bdSLois Curfman McInnes   } else {
849f830108cSBarry Smith     PetscOps       *Abops;
850f830108cSBarry Smith     struct _MatOps *Aops;
851f830108cSBarry Smith 
8523501a2bdSLois Curfman McInnes     /* This isn't really an in-place transpose, but free data struct from a */
8530452661fSBarry Smith     PetscFree(a->rowners);
8543501a2bdSLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
8553501a2bdSLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
8563501a2bdSLois Curfman McInnes     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
8570452661fSBarry Smith     PetscFree(a);
858f830108cSBarry Smith 
859f830108cSBarry Smith     /*
860f830108cSBarry Smith          This is horrible, horrible code. We need to keep the
861f830108cSBarry Smith       A pointers for the bops and ops but copy everything
862f830108cSBarry Smith       else from C.
863f830108cSBarry Smith     */
864f830108cSBarry Smith     Abops = A->bops;
865f830108cSBarry Smith     Aops  = A->ops;
866f09e8eb9SSatish Balay     PetscMemcpy(A,B,sizeof(struct _p_Mat));
867f830108cSBarry Smith     A->bops = Abops;
868f830108cSBarry Smith     A->ops  = Aops;
869f830108cSBarry Smith 
8700452661fSBarry Smith     PetscHeaderDestroy(B);
8713501a2bdSLois Curfman McInnes   }
8723a40ed3dSBarry Smith   PetscFunctionReturn(0);
873096963f5SLois Curfman McInnes }
874096963f5SLois Curfman McInnes 
875eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
8765615d1e5SSatish Balay #undef __FUNC__
8775615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense"
8788f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA)
87944cd7ae7SLois Curfman McInnes {
88044cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense *) inA->data;
88144cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->A->data;
88244cd7ae7SLois Curfman McInnes   int          one = 1, nz;
88344cd7ae7SLois Curfman McInnes 
8843a40ed3dSBarry Smith   PetscFunctionBegin;
88544cd7ae7SLois Curfman McInnes   nz = a->m*a->n;
88644cd7ae7SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
88744cd7ae7SLois Curfman McInnes   PLogFlops(nz);
8883a40ed3dSBarry Smith   PetscFunctionReturn(0);
88944cd7ae7SLois Curfman McInnes }
89044cd7ae7SLois Curfman McInnes 
8913d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat,Mat *,int);
8928965ea79SLois Curfman McInnes 
8938965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
89439ddd567SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIDense,
89539ddd567SLois Curfman McInnes        MatGetRow_MPIDense,MatRestoreRow_MPIDense,
89639ddd567SLois Curfman McInnes        MatMult_MPIDense,MatMultAdd_MPIDense,
897096963f5SLois Curfman McInnes        MatMultTrans_MPIDense,MatMultTransAdd_MPIDense,
898e7ca0642SLois Curfman McInnes        0,0,
89939ddd567SLois Curfman McInnes        0,0,
9008c469469SLois Curfman McInnes        0,0,
9018aaee692SLois Curfman McInnes        0,MatTranspose_MPIDense,
90239ddd567SLois Curfman McInnes        MatGetInfo_MPIDense,0,
903096963f5SLois Curfman McInnes        MatGetDiagonal_MPIDense,0,MatNorm_MPIDense,
90439ddd567SLois Curfman McInnes        MatAssemblyBegin_MPIDense,MatAssemblyEnd_MPIDense,
9058965ea79SLois Curfman McInnes        0,
90639ddd567SLois Curfman McInnes        MatSetOption_MPIDense,MatZeroEntries_MPIDense,MatZeroRows_MPIDense,
9073b2fbd54SBarry Smith        0,0,
9088aaee692SLois Curfman McInnes        0,0,
90939ddd567SLois Curfman McInnes        MatGetSize_MPIDense,MatGetLocalSize_MPIDense,
91039ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
911ff14e315SSatish Balay        0,0, MatGetArray_MPIDense, MatRestoreArray_MPIDense,
91294a9d846SBarry Smith        MatConvertSameType_MPIDense,
913b49de8d1SLois Curfman McInnes        0,0,0,0,0,
9144e220ebcSLois Curfman McInnes        0,0,MatGetValues_MPIDense,0,0,MatScale_MPIDense,
9154e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_MPIDense};
9168965ea79SLois Curfman McInnes 
9175615d1e5SSatish Balay #undef __FUNC__
9185615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense"
9198965ea79SLois Curfman McInnes /*@C
92039ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
9218965ea79SLois Curfman McInnes 
922db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
923db81eaa0SLois Curfman McInnes 
9248965ea79SLois Curfman McInnes    Input Parameters:
925db81eaa0SLois Curfman McInnes +  comm - MPI communicator
9268965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
927db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
9288965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
929db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
930db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
931dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
9328965ea79SLois Curfman McInnes 
9338965ea79SLois Curfman McInnes    Output Parameter:
934477f1c0bSLois Curfman McInnes .  A - the matrix
9358965ea79SLois Curfman McInnes 
936b259b22eSLois Curfman McInnes    Notes:
93739ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
93839ddd567SLois Curfman McInnes    storage by columns.
9398965ea79SLois Curfman McInnes 
94018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
94118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
942b4fd4287SBarry Smith    set data=PETSC_NULL.
94318f449edSLois Curfman McInnes 
9448965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
9458965ea79SLois Curfman McInnes    (possibly both).
9468965ea79SLois Curfman McInnes 
9473501a2bdSLois Curfman McInnes    Currently, the only parallel dense matrix decomposition is by rows,
9483501a2bdSLois Curfman McInnes    so that n=N and each submatrix owns all of the global columns.
9493501a2bdSLois Curfman McInnes 
95039ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel
9518965ea79SLois Curfman McInnes 
95239ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
9538965ea79SLois Curfman McInnes @*/
954477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
9558965ea79SLois Curfman McInnes {
9568965ea79SLois Curfman McInnes   Mat          mat;
95739ddd567SLois Curfman McInnes   Mat_MPIDense *a;
95825cdf11fSBarry Smith   int          ierr, i,flg;
9598965ea79SLois Curfman McInnes 
9603a40ed3dSBarry Smith   PetscFunctionBegin;
961ed2daf61SLois Curfman McInnes   /* Note:  For now, when data is specified above, this assumes the user correctly
962ed2daf61SLois Curfman McInnes    allocates the local dense storage space.  We should add error checking. */
96318f449edSLois Curfman McInnes 
964477f1c0bSLois Curfman McInnes   *A = 0;
965f830108cSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,comm,MatDestroy,MatView);
9668965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
9670452661fSBarry Smith   mat->data       = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
968f830108cSBarry Smith   PetscMemcpy(mat->ops,&MatOps,sizeof(struct _MatOps));
969e1311b90SBarry Smith   mat->ops->destroy    = MatDestroy_MPIDense;
970e1311b90SBarry Smith   mat->ops->view       = MatView_MPIDense;
9718965ea79SLois Curfman McInnes   mat->factor          = 0;
97290f02eecSBarry Smith   mat->mapping         = 0;
9738965ea79SLois Curfman McInnes 
974622d7880SLois Curfman McInnes   a->factor       = 0;
975e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
9768965ea79SLois Curfman McInnes   MPI_Comm_rank(comm,&a->rank);
9778965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&a->size);
9788965ea79SLois Curfman McInnes 
979ca161407SBarry Smith   if (M == PETSC_DECIDE) {ierr = MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);}
9808965ea79SLois Curfman McInnes   if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);}
98139ddd567SLois Curfman McInnes 
98239ddd567SLois Curfman McInnes   /* each row stores all columns */
98339ddd567SLois Curfman McInnes   if (N == PETSC_DECIDE) N = n;
984d7e8b826SBarry Smith   if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);}
985a8c6a408SBarry Smith   /*  if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */
986aca0ad90SLois Curfman McInnes   a->N = mat->N = N;
987aca0ad90SLois Curfman McInnes   a->M = mat->M = M;
988aca0ad90SLois Curfman McInnes   a->m = mat->m = m;
989aca0ad90SLois Curfman McInnes   a->n = mat->n = n;
9908965ea79SLois Curfman McInnes 
9918965ea79SLois Curfman McInnes   /* build local table of row and column ownerships */
992d7e8b826SBarry Smith   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners);
993d7e8b826SBarry Smith   a->cowners = a->rowners + a->size + 1;
994f09e8eb9SSatish Balay   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
995ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
9968965ea79SLois Curfman McInnes   a->rowners[0] = 0;
9978965ea79SLois Curfman McInnes   for ( i=2; i<=a->size; i++ ) {
9988965ea79SLois Curfman McInnes     a->rowners[i] += a->rowners[i-1];
9998965ea79SLois Curfman McInnes   }
10008965ea79SLois Curfman McInnes   a->rstart = a->rowners[a->rank];
10018965ea79SLois Curfman McInnes   a->rend   = a->rowners[a->rank+1];
1002ca161407SBarry Smith   ierr      = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
1003d7e8b826SBarry Smith   a->cowners[0] = 0;
1004d7e8b826SBarry Smith   for ( i=2; i<=a->size; i++ ) {
1005d7e8b826SBarry Smith     a->cowners[i] += a->cowners[i-1];
1006d7e8b826SBarry Smith   }
10078965ea79SLois Curfman McInnes 
1008029af93fSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A); CHKERRQ(ierr);
10098965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
10108965ea79SLois Curfman McInnes 
10118965ea79SLois Curfman McInnes   /* build cache for off array entries formed */
10128965ea79SLois Curfman McInnes   ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr);
10138965ea79SLois Curfman McInnes 
10148965ea79SLois Curfman McInnes   /* stuff used for matrix vector multiply */
10158965ea79SLois Curfman McInnes   a->lvec        = 0;
10168965ea79SLois Curfman McInnes   a->Mvctx       = 0;
101739b7565bSBarry Smith   a->roworiented = 1;
10188965ea79SLois Curfman McInnes 
1019477f1c0bSLois Curfman McInnes   *A = mat;
102025cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
102125cdf11fSBarry Smith   if (flg) {
10228c469469SLois Curfman McInnes     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
10238c469469SLois Curfman McInnes   }
10243a40ed3dSBarry Smith   PetscFunctionReturn(0);
10258965ea79SLois Curfman McInnes }
10268965ea79SLois Curfman McInnes 
10275615d1e5SSatish Balay #undef __FUNC__
10285615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIDense"
10293d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat A,Mat *newmat,int cpvalues)
10308965ea79SLois Curfman McInnes {
10318965ea79SLois Curfman McInnes   Mat          mat;
10323501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data;
103339ddd567SLois Curfman McInnes   int          ierr;
10342ba99913SLois Curfman McInnes   FactorCtx    *factor;
10358965ea79SLois Curfman McInnes 
10363a40ed3dSBarry Smith   PetscFunctionBegin;
10378965ea79SLois Curfman McInnes   *newmat       = 0;
1038f830108cSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,A->comm,MatDestroy,MatView);
10398965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
10400452661fSBarry Smith   mat->data      = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
1041f830108cSBarry Smith   PetscMemcpy(mat->ops,&MatOps,sizeof(struct _MatOps));
1042e1311b90SBarry Smith   mat->ops->destroy   = MatDestroy_MPIDense;
1043e1311b90SBarry Smith   mat->ops->view      = MatView_MPIDense;
10443501a2bdSLois Curfman McInnes   mat->factor         = A->factor;
1045c456f294SBarry Smith   mat->assembled      = PETSC_TRUE;
10468965ea79SLois Curfman McInnes 
104744cd7ae7SLois Curfman McInnes   a->m = mat->m = oldmat->m;
104844cd7ae7SLois Curfman McInnes   a->n = mat->n = oldmat->n;
104944cd7ae7SLois Curfman McInnes   a->M = mat->M = oldmat->M;
105044cd7ae7SLois Curfman McInnes   a->N = mat->N = oldmat->N;
10512ba99913SLois Curfman McInnes   if (oldmat->factor) {
10522ba99913SLois Curfman McInnes     a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx)); CHKPTRQ(factor);
10532ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
10542ba99913SLois Curfman McInnes   } else a->factor = 0;
10558965ea79SLois Curfman McInnes 
10568965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
10578965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
10588965ea79SLois Curfman McInnes   a->size         = oldmat->size;
10598965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1060e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
10618965ea79SLois Curfman McInnes 
10620452661fSBarry Smith   a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners);
1063f09e8eb9SSatish Balay   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
10648965ea79SLois Curfman McInnes   PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));
10658965ea79SLois Curfman McInnes   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
10668965ea79SLois Curfman McInnes 
10678965ea79SLois Curfman McInnes   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
10688965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->lvec);
106955659b69SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
10708965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->Mvctx);
10718965ea79SLois Curfman McInnes   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
10728965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
10738965ea79SLois Curfman McInnes   *newmat = mat;
10743a40ed3dSBarry Smith   PetscFunctionReturn(0);
10758965ea79SLois Curfman McInnes }
10768965ea79SLois Curfman McInnes 
107777c4ece6SBarry Smith #include "sys.h"
10788965ea79SLois Curfman McInnes 
10795615d1e5SSatish Balay #undef __FUNC__
10805615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile"
108190ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat)
108290ace30eSBarry Smith {
108340011551SBarry Smith   int        *rowners, i,size,rank,m,ierr,nz,j;
108490ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
108590ace30eSBarry Smith   MPI_Status status;
108690ace30eSBarry Smith 
10873a40ed3dSBarry Smith   PetscFunctionBegin;
108890ace30eSBarry Smith   MPI_Comm_rank(comm,&rank);
108990ace30eSBarry Smith   MPI_Comm_size(comm,&size);
109090ace30eSBarry Smith 
109190ace30eSBarry Smith   /* determine ownership of all rows */
109290ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
109390ace30eSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
1094ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
109590ace30eSBarry Smith   rowners[0] = 0;
109690ace30eSBarry Smith   for ( i=2; i<=size; i++ ) {
109790ace30eSBarry Smith     rowners[i] += rowners[i-1];
109890ace30eSBarry Smith   }
109990ace30eSBarry Smith 
110090ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
110190ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array); CHKERRQ(ierr);
110290ace30eSBarry Smith 
110390ace30eSBarry Smith   if (!rank) {
110490ace30eSBarry Smith     vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
110590ace30eSBarry Smith 
110690ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
11070752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR); CHKERRQ(ierr);
110890ace30eSBarry Smith 
110990ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
111090ace30eSBarry Smith     vals_ptr = vals;
111190ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
111290ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
111390ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
111490ace30eSBarry Smith       }
111590ace30eSBarry Smith     }
111690ace30eSBarry Smith 
111790ace30eSBarry Smith     /* read in other processors and ship out */
111890ace30eSBarry Smith     for ( i=1; i<size; i++ ) {
111990ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
11200752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
1121ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
112290ace30eSBarry Smith     }
11233a40ed3dSBarry Smith   } else {
112490ace30eSBarry Smith     /* receive numeric values */
112590ace30eSBarry Smith     vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
112690ace30eSBarry Smith 
112790ace30eSBarry Smith     /* receive message of values*/
1128ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
112990ace30eSBarry Smith 
113090ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
113190ace30eSBarry Smith     vals_ptr = vals;
113290ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
113390ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
113490ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
113590ace30eSBarry Smith       }
113690ace30eSBarry Smith     }
113790ace30eSBarry Smith   }
113890ace30eSBarry Smith   PetscFree(rowners);
113990ace30eSBarry Smith   PetscFree(vals);
11406d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11416d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11423a40ed3dSBarry Smith   PetscFunctionReturn(0);
114390ace30eSBarry Smith }
114490ace30eSBarry Smith 
114590ace30eSBarry Smith 
11465615d1e5SSatish Balay #undef __FUNC__
11475615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense"
114819bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
11498965ea79SLois Curfman McInnes {
11508965ea79SLois Curfman McInnes   Mat          A;
11518965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
115219bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
11538965ea79SLois Curfman McInnes   MPI_Status   status;
11548965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
11558965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
115619bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
11573a40ed3dSBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
11588965ea79SLois Curfman McInnes 
11593a40ed3dSBarry Smith   PetscFunctionBegin;
11608965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank);
11618965ea79SLois Curfman McInnes   if (!rank) {
116290ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
11630752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr);
1164a8c6a408SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object");
11658965ea79SLois Curfman McInnes   }
11668965ea79SLois Curfman McInnes 
1167ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
116890ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
116990ace30eSBarry Smith 
117090ace30eSBarry Smith   /*
117190ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
117290ace30eSBarry Smith   */
117390ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
11743a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
11753a40ed3dSBarry Smith     PetscFunctionReturn(0);
117690ace30eSBarry Smith   }
117790ace30eSBarry Smith 
11788965ea79SLois Curfman McInnes   /* determine ownership of all rows */
11798965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
11800452661fSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
1181ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
11828965ea79SLois Curfman McInnes   rowners[0] = 0;
11838965ea79SLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
11848965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
11858965ea79SLois Curfman McInnes   }
11868965ea79SLois Curfman McInnes   rstart = rowners[rank];
11878965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
11888965ea79SLois Curfman McInnes 
11898965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
11900452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
11918965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
11928965ea79SLois Curfman McInnes   if (!rank) {
11930452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
11940752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
11950452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts);
11968965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1197ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
11980452661fSBarry Smith     PetscFree(sndcounts);
1199ca161407SBarry Smith   } else {
1200ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr);
12018965ea79SLois Curfman McInnes   }
12028965ea79SLois Curfman McInnes 
12038965ea79SLois Curfman McInnes   if (!rank) {
12048965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
12050452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz);
1206cddf8d76SBarry Smith     PetscMemzero(procsnz,size*sizeof(int));
12078965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
12088965ea79SLois Curfman McInnes       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
12098965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
12108965ea79SLois Curfman McInnes       }
12118965ea79SLois Curfman McInnes     }
12120452661fSBarry Smith     PetscFree(rowlengths);
12138965ea79SLois Curfman McInnes 
12148965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
12158965ea79SLois Curfman McInnes     maxnz = 0;
12168965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
12170452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
12188965ea79SLois Curfman McInnes     }
12190452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols);
12208965ea79SLois Curfman McInnes 
12218965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
12228965ea79SLois Curfman McInnes     nz = procsnz[0];
12230452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
12240752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr);
12258965ea79SLois Curfman McInnes 
12268965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
12278965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
12288965ea79SLois Curfman McInnes       nz   = procsnz[i];
12290752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
1230ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
12318965ea79SLois Curfman McInnes     }
12320452661fSBarry Smith     PetscFree(cols);
12333a40ed3dSBarry Smith   } else {
12348965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
12358965ea79SLois Curfman McInnes     nz = 0;
12368965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
12378965ea79SLois Curfman McInnes       nz += ourlens[i];
12388965ea79SLois Curfman McInnes     }
12390452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
12408965ea79SLois Curfman McInnes 
12418965ea79SLois Curfman McInnes     /* receive message of column indices*/
1242ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1243ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
1244a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
12458965ea79SLois Curfman McInnes   }
12468965ea79SLois Curfman McInnes 
12478965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1248cddf8d76SBarry Smith   PetscMemzero(offlens,m*sizeof(int));
12498965ea79SLois Curfman McInnes   jj = 0;
12508965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
12518965ea79SLois Curfman McInnes     for ( j=0; j<ourlens[i]; j++ ) {
12528965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
12538965ea79SLois Curfman McInnes       jj++;
12548965ea79SLois Curfman McInnes     }
12558965ea79SLois Curfman McInnes   }
12568965ea79SLois Curfman McInnes 
12578965ea79SLois Curfman McInnes   /* create our matrix */
12588965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
12598965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
12608965ea79SLois Curfman McInnes   }
1261b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
12628965ea79SLois Curfman McInnes   A = *newmat;
12638965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
12648965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
12658965ea79SLois Curfman McInnes   }
12668965ea79SLois Curfman McInnes 
12678965ea79SLois Curfman McInnes   if (!rank) {
12680452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
12698965ea79SLois Curfman McInnes 
12708965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
12718965ea79SLois Curfman McInnes     nz = procsnz[0];
12720752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
12738965ea79SLois Curfman McInnes 
12748965ea79SLois Curfman McInnes     /* insert into matrix */
12758965ea79SLois Curfman McInnes     jj      = rstart;
12768965ea79SLois Curfman McInnes     smycols = mycols;
12778965ea79SLois Curfman McInnes     svals   = vals;
12788965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
12798965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
12808965ea79SLois Curfman McInnes       smycols += ourlens[i];
12818965ea79SLois Curfman McInnes       svals   += ourlens[i];
12828965ea79SLois Curfman McInnes       jj++;
12838965ea79SLois Curfman McInnes     }
12848965ea79SLois Curfman McInnes 
12858965ea79SLois Curfman McInnes     /* read in other processors and ship out */
12868965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
12878965ea79SLois Curfman McInnes       nz   = procsnz[i];
12880752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
1289ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
12908965ea79SLois Curfman McInnes     }
12910452661fSBarry Smith     PetscFree(procsnz);
12923a40ed3dSBarry Smith   } else {
12938965ea79SLois Curfman McInnes     /* receive numeric values */
12940452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
12958965ea79SLois Curfman McInnes 
12968965ea79SLois Curfman McInnes     /* receive message of values*/
1297ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1298ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
1299a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
13008965ea79SLois Curfman McInnes 
13018965ea79SLois Curfman McInnes     /* insert into matrix */
13028965ea79SLois Curfman McInnes     jj      = rstart;
13038965ea79SLois Curfman McInnes     smycols = mycols;
13048965ea79SLois Curfman McInnes     svals   = vals;
13058965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
13068965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13078965ea79SLois Curfman McInnes       smycols += ourlens[i];
13088965ea79SLois Curfman McInnes       svals   += ourlens[i];
13098965ea79SLois Curfman McInnes       jj++;
13108965ea79SLois Curfman McInnes     }
13118965ea79SLois Curfman McInnes   }
13120452661fSBarry Smith   PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners);
13138965ea79SLois Curfman McInnes 
13146d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13156d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13163a40ed3dSBarry Smith   PetscFunctionReturn(0);
13178965ea79SLois Curfman McInnes }
131890ace30eSBarry Smith 
131990ace30eSBarry Smith 
132090ace30eSBarry Smith 
132190ace30eSBarry Smith 
132290ace30eSBarry Smith 
1323