xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 61b13de0c92b07d489559c7f7ccd03b9191851f1)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*61b13de0SBarry Smith static char vcid[] = "$Id: mpidense.c,v 1.93 1998/07/13 20:34:24 bsmith Exp bsmith $";
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);
501*61b13de0SBarry Smith   if (mat->rmap) {
502*61b13de0SBarry Smith     ierr = MapDestroy(mat->rmap);CHKERRQ(ierr);
503*61b13de0SBarry Smith   }
504*61b13de0SBarry Smith   if (mat->cmap) {
505*61b13de0SBarry Smith     ierr = MapDestroy(mat->cmap);CHKERRQ(ierr);
50690f02eecSBarry Smith   }
5078965ea79SLois Curfman McInnes   PLogObjectDestroy(mat);
5080452661fSBarry Smith   PetscHeaderDestroy(mat);
5093a40ed3dSBarry Smith   PetscFunctionReturn(0);
5108965ea79SLois Curfman McInnes }
51139ddd567SLois Curfman McInnes 
5128965ea79SLois Curfman McInnes #include "pinclude/pviewer.h"
5138965ea79SLois Curfman McInnes 
5145615d1e5SSatish Balay #undef __FUNC__
5155615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary"
51639ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
5178965ea79SLois Curfman McInnes {
51839ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
5198965ea79SLois Curfman McInnes   int          ierr;
5207056b6fcSBarry Smith 
5213a40ed3dSBarry Smith   PetscFunctionBegin;
52239ddd567SLois Curfman McInnes   if (mdn->size == 1) {
52339ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5248965ea79SLois Curfman McInnes   }
525a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported");
5263a40ed3dSBarry Smith   PetscFunctionReturn(0);
5278965ea79SLois Curfman McInnes }
5288965ea79SLois Curfman McInnes 
5295615d1e5SSatish Balay #undef __FUNC__
5305615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_ASCII"
53139ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer)
5328965ea79SLois Curfman McInnes {
53339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
53439ddd567SLois Curfman McInnes   int          ierr, format;
5358965ea79SLois Curfman McInnes   FILE         *fd;
53619bcc07fSBarry Smith   ViewerType   vtype;
5378965ea79SLois Curfman McInnes 
5383a40ed3dSBarry Smith   PetscFunctionBegin;
5393a40ed3dSBarry Smith   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
54090ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
54190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
542639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
5434e220ebcSLois Curfman McInnes     int rank;
5444e220ebcSLois Curfman McInnes     MatInfo info;
545096963f5SLois Curfman McInnes     MPI_Comm_rank(mat->comm,&rank);
5464e220ebcSLois Curfman McInnes     ierr = MatGetInfo(mat,MAT_LOCAL,&info);
54777c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
5484e220ebcSLois Curfman McInnes       fprintf(fd,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m,
5494e220ebcSLois Curfman McInnes          (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
550096963f5SLois Curfman McInnes       fflush(fd);
55177c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
5523501a2bdSLois Curfman McInnes     ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr);
5533a40ed3dSBarry Smith     PetscFunctionReturn(0);
5543501a2bdSLois Curfman McInnes   }
55502cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_INFO) {
5563a40ed3dSBarry Smith     PetscFunctionReturn(0);
5578965ea79SLois Curfman McInnes   }
55819bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
55977c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
56039ddd567SLois Curfman McInnes     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n",
56139ddd567SLois Curfman McInnes              mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n);
56239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5638965ea79SLois Curfman McInnes     fflush(fd);
56477c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
5653a40ed3dSBarry Smith   } else {
56639ddd567SLois Curfman McInnes     int size = mdn->size, rank = mdn->rank;
5678965ea79SLois Curfman McInnes     if (size == 1) {
56839ddd567SLois Curfman McInnes       ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5693a40ed3dSBarry Smith     } else {
5708965ea79SLois Curfman McInnes       /* assemble the entire matrix onto first processor. */
5718965ea79SLois Curfman McInnes       Mat          A;
57239ddd567SLois Curfman McInnes       int          M = mdn->M, N = mdn->N,m,row,i, nz, *cols;
57339ddd567SLois Curfman McInnes       Scalar       *vals;
57439ddd567SLois Curfman McInnes       Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data;
5758965ea79SLois Curfman McInnes 
5768965ea79SLois Curfman McInnes       if (!rank) {
5770513a670SBarry Smith         ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr);
5783a40ed3dSBarry Smith       } else {
5790513a670SBarry Smith         ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr);
5808965ea79SLois Curfman McInnes       }
5818965ea79SLois Curfman McInnes       PLogObjectParent(mat,A);
5828965ea79SLois Curfman McInnes 
58339ddd567SLois Curfman McInnes       /* Copy the matrix ... This isn't the most efficient means,
58439ddd567SLois Curfman McInnes          but it's quick for now */
58539ddd567SLois Curfman McInnes       row = mdn->rstart; m = Amdn->m;
5868965ea79SLois Curfman McInnes       for ( i=0; i<m; i++ ) {
58739ddd567SLois Curfman McInnes         ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
58839ddd567SLois Curfman McInnes         ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr);
58939ddd567SLois Curfman McInnes         ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
59039ddd567SLois Curfman McInnes         row++;
5918965ea79SLois Curfman McInnes       }
5928965ea79SLois Curfman McInnes 
5936d4a8577SBarry Smith       ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5946d4a8577SBarry Smith       ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5958965ea79SLois Curfman McInnes       if (!rank) {
59639ddd567SLois Curfman McInnes         ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr);
5978965ea79SLois Curfman McInnes       }
5988965ea79SLois Curfman McInnes       ierr = MatDestroy(A); CHKERRQ(ierr);
5998965ea79SLois Curfman McInnes     }
6008965ea79SLois Curfman McInnes   }
6013a40ed3dSBarry Smith   PetscFunctionReturn(0);
6028965ea79SLois Curfman McInnes }
6038965ea79SLois Curfman McInnes 
6045615d1e5SSatish Balay #undef __FUNC__
6055615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense"
606e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer)
6078965ea79SLois Curfman McInnes {
60839ddd567SLois Curfman McInnes   int          ierr;
609bcd2baecSBarry Smith   ViewerType   vtype;
6108965ea79SLois Curfman McInnes 
611bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
612bcd2baecSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
61339ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
6143a40ed3dSBarry Smith   } else if (vtype == ASCII_FILES_VIEWER) {
61539ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
6163a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
6173a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6185cd90555SBarry Smith   } else {
6195cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
6208965ea79SLois Curfman McInnes   }
6213a40ed3dSBarry Smith   PetscFunctionReturn(0);
6228965ea79SLois Curfman McInnes }
6238965ea79SLois Curfman McInnes 
6245615d1e5SSatish Balay #undef __FUNC__
6255615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense"
6268f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6278965ea79SLois Curfman McInnes {
6283501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6293501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6304e220ebcSLois Curfman McInnes   int          ierr;
6314e220ebcSLois Curfman McInnes   double       isend[5], irecv[5];
6328965ea79SLois Curfman McInnes 
6333a40ed3dSBarry Smith   PetscFunctionBegin;
6344e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
6354e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
6364e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
6374e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
6384e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6394e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info); CHKERRQ(ierr);
6404e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6414e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6428965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6434e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6444e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6454e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6464e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6474e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6488965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
649f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
6504e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6514e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6524e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6534e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6544e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6558965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
656f7cdd7c9SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
6574e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6584e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6594e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6604e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6614e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6628965ea79SLois Curfman McInnes   }
6634e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6644e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6654e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6663a40ed3dSBarry Smith   PetscFunctionReturn(0);
6678965ea79SLois Curfman McInnes }
6688965ea79SLois Curfman McInnes 
6698c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*);
6708aaee692SLois Curfman McInnes    extern int MatLUFactorNumeric_MPIDense(Mat,Mat*);
6718aaee692SLois Curfman McInnes    extern int MatLUFactor_MPIDense(Mat,IS,IS,double);
6728aaee692SLois Curfman McInnes    extern int MatSolve_MPIDense(Mat,Vec,Vec);
6738c469469SLois Curfman McInnes    extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec);
6748aaee692SLois Curfman McInnes    extern int MatSolveTrans_MPIDense(Mat,Vec,Vec);
6758aaee692SLois Curfman McInnes    extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */
6768aaee692SLois Curfman McInnes 
6775615d1e5SSatish Balay #undef __FUNC__
6785615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense"
6798f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6808965ea79SLois Curfman McInnes {
68139ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
6828965ea79SLois Curfman McInnes 
6833a40ed3dSBarry Smith   PetscFunctionBegin;
6846d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6856d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
686c2653b3dSLois Curfman McInnes       op == MAT_NEW_NONZERO_LOCATION_ERROR ||
68796854ed6SLois Curfman McInnes       op == MAT_NEW_NONZERO_ALLOCATION_ERROR ||
688219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
689219d9a1aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED) {
690b1fbbac0SLois Curfman McInnes         MatSetOption(a->A,op);
691b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
692aeafbbfcSLois Curfman McInnes         a->roworiented = 1;
6938965ea79SLois Curfman McInnes         MatSetOption(a->A,op);
694b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
695219d9a1aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
6966d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
6976d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
698b51ba29fSSatish Balay              op == MAT_YES_NEW_DIAGONALS ||
699b51ba29fSSatish Balay              op == MAT_USE_HASH_TABLE) {
700981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
7013a40ed3dSBarry Smith   } else if (op == MAT_COLUMN_ORIENTED) {
7023a40ed3dSBarry Smith     a->roworiented = 0; MatSetOption(a->A,op);
7033a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS) {
7043a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7053a40ed3dSBarry Smith   } else {
7063a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7073a40ed3dSBarry Smith   }
7083a40ed3dSBarry Smith   PetscFunctionReturn(0);
7098965ea79SLois Curfman McInnes }
7108965ea79SLois Curfman McInnes 
7115615d1e5SSatish Balay #undef __FUNC__
7125615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense"
7138f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n)
7148965ea79SLois Curfman McInnes {
7153501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7163a40ed3dSBarry Smith 
7173a40ed3dSBarry Smith   PetscFunctionBegin;
7188965ea79SLois Curfman McInnes   *m = mat->M; *n = mat->N;
7193a40ed3dSBarry Smith   PetscFunctionReturn(0);
7208965ea79SLois Curfman McInnes }
7218965ea79SLois Curfman McInnes 
7225615d1e5SSatish Balay #undef __FUNC__
7235615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense"
7248f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n)
7258965ea79SLois Curfman McInnes {
7263501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7273a40ed3dSBarry Smith 
7283a40ed3dSBarry Smith   PetscFunctionBegin;
7298965ea79SLois Curfman McInnes   *m = mat->m; *n = mat->N;
7303a40ed3dSBarry Smith   PetscFunctionReturn(0);
7318965ea79SLois Curfman McInnes }
7328965ea79SLois Curfman McInnes 
7335615d1e5SSatish Balay #undef __FUNC__
7345615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense"
7358f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
7368965ea79SLois Curfman McInnes {
7373501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7383a40ed3dSBarry Smith 
7393a40ed3dSBarry Smith   PetscFunctionBegin;
7408965ea79SLois Curfman McInnes   *m = mat->rstart; *n = mat->rend;
7413a40ed3dSBarry Smith   PetscFunctionReturn(0);
7428965ea79SLois Curfman McInnes }
7438965ea79SLois Curfman McInnes 
7445615d1e5SSatish Balay #undef __FUNC__
7455615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense"
7468f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
7478965ea79SLois Curfman McInnes {
7483501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
7493a40ed3dSBarry Smith   int          lrow, rstart = mat->rstart, rend = mat->rend,ierr;
7508965ea79SLois Curfman McInnes 
7513a40ed3dSBarry Smith   PetscFunctionBegin;
752a8c6a408SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows")
7538965ea79SLois Curfman McInnes   lrow = row - rstart;
7543a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7553a40ed3dSBarry Smith   PetscFunctionReturn(0);
7568965ea79SLois Curfman McInnes }
7578965ea79SLois Curfman McInnes 
7585615d1e5SSatish Balay #undef __FUNC__
7595615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense"
7608f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
7618965ea79SLois Curfman McInnes {
7623a40ed3dSBarry Smith   PetscFunctionBegin;
7630452661fSBarry Smith   if (idx) PetscFree(*idx);
7640452661fSBarry Smith   if (v) PetscFree(*v);
7653a40ed3dSBarry Smith   PetscFunctionReturn(0);
7668965ea79SLois Curfman McInnes }
7678965ea79SLois Curfman McInnes 
7685615d1e5SSatish Balay #undef __FUNC__
7695615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense"
7708f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm)
771096963f5SLois Curfman McInnes {
7723501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
7733501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
7743501a2bdSLois Curfman McInnes   int          ierr, i, j;
7753501a2bdSLois Curfman McInnes   double       sum = 0.0;
7763501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
7773501a2bdSLois Curfman McInnes 
7783a40ed3dSBarry Smith   PetscFunctionBegin;
7793501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
7803501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm); CHKERRQ(ierr);
7813501a2bdSLois Curfman McInnes   } else {
7823501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
7833501a2bdSLois Curfman McInnes       for (i=0; i<mat->n*mat->m; i++ ) {
7843a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
785e20fef11SSatish Balay         sum += PetscReal(PetscConj(*v)*(*v)); v++;
7863501a2bdSLois Curfman McInnes #else
7873501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7883501a2bdSLois Curfman McInnes #endif
7893501a2bdSLois Curfman McInnes       }
790ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
7913501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
7923501a2bdSLois Curfman McInnes       PLogFlops(2*mat->n*mat->m);
7933a40ed3dSBarry Smith     } else if (type == NORM_1) {
7943501a2bdSLois Curfman McInnes       double *tmp, *tmp2;
7950452661fSBarry Smith       tmp  = (double *) PetscMalloc( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp);
7963501a2bdSLois Curfman McInnes       tmp2 = tmp + mdn->N;
797cddf8d76SBarry Smith       PetscMemzero(tmp,2*mdn->N*sizeof(double));
798096963f5SLois Curfman McInnes       *norm = 0.0;
7993501a2bdSLois Curfman McInnes       v = mat->v;
8003501a2bdSLois Curfman McInnes       for ( j=0; j<mat->n; j++ ) {
8013501a2bdSLois Curfman McInnes         for ( i=0; i<mat->m; i++ ) {
80267e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8033501a2bdSLois Curfman McInnes         }
8043501a2bdSLois Curfman McInnes       }
805ca161407SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr);
8063501a2bdSLois Curfman McInnes       for ( j=0; j<mdn->N; j++ ) {
8073501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
8083501a2bdSLois Curfman McInnes       }
8090452661fSBarry Smith       PetscFree(tmp);
8103501a2bdSLois Curfman McInnes       PLogFlops(mat->n*mat->m);
8113a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
8123501a2bdSLois Curfman McInnes       double ntemp;
8133501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr);
814ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr);
8153a40ed3dSBarry Smith     } else {
816a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
8173501a2bdSLois Curfman McInnes     }
8183501a2bdSLois Curfman McInnes   }
8193a40ed3dSBarry Smith   PetscFunctionReturn(0);
8203501a2bdSLois Curfman McInnes }
8213501a2bdSLois Curfman McInnes 
8225615d1e5SSatish Balay #undef __FUNC__
8235615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense"
8248f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8253501a2bdSLois Curfman McInnes {
8263501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
8273501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data;
8283501a2bdSLois Curfman McInnes   Mat          B;
8293501a2bdSLois Curfman McInnes   int          M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart;
8303501a2bdSLois Curfman McInnes   int          j, i, ierr;
8313501a2bdSLois Curfman McInnes   Scalar       *v;
8323501a2bdSLois Curfman McInnes 
8333a40ed3dSBarry Smith   PetscFunctionBegin;
8347056b6fcSBarry Smith   if (matout == PETSC_NULL && M != N) {
835a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place");
8367056b6fcSBarry Smith   }
8377056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8383501a2bdSLois Curfman McInnes 
8393501a2bdSLois Curfman McInnes   m = Aloc->m; n = Aloc->n; v = Aloc->v;
8400452661fSBarry Smith   rwork = (int *) PetscMalloc(n*sizeof(int)); CHKPTRQ(rwork);
8413501a2bdSLois Curfman McInnes   for ( j=0; j<n; j++ ) {
8423501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8433501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr);
8443501a2bdSLois Curfman McInnes     v   += m;
8453501a2bdSLois Curfman McInnes   }
8460452661fSBarry Smith   PetscFree(rwork);
8476d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8486d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8493638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
8503501a2bdSLois Curfman McInnes     *matout = B;
8513501a2bdSLois Curfman McInnes   } else {
852f830108cSBarry Smith     PetscOps *Abops;
85309dc0095SBarry Smith     MatOps   Aops;
854f830108cSBarry Smith 
8553501a2bdSLois Curfman McInnes     /* This isn't really an in-place transpose, but free data struct from a */
8560452661fSBarry Smith     PetscFree(a->rowners);
8573501a2bdSLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
8583501a2bdSLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
8593501a2bdSLois Curfman McInnes     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
8600452661fSBarry Smith     PetscFree(a);
861f830108cSBarry Smith 
862f830108cSBarry Smith     /*
863f830108cSBarry Smith          This is horrible, horrible code. We need to keep the
864f830108cSBarry Smith       A pointers for the bops and ops but copy everything
865f830108cSBarry Smith       else from C.
866f830108cSBarry Smith     */
867f830108cSBarry Smith     Abops = A->bops;
868f830108cSBarry Smith     Aops  = A->ops;
869f09e8eb9SSatish Balay     PetscMemcpy(A,B,sizeof(struct _p_Mat));
870f830108cSBarry Smith     A->bops = Abops;
871f830108cSBarry Smith     A->ops  = Aops;
872f830108cSBarry Smith 
8730452661fSBarry Smith     PetscHeaderDestroy(B);
8743501a2bdSLois Curfman McInnes   }
8753a40ed3dSBarry Smith   PetscFunctionReturn(0);
876096963f5SLois Curfman McInnes }
877096963f5SLois Curfman McInnes 
878eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
8795615d1e5SSatish Balay #undef __FUNC__
8805615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense"
8818f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA)
88244cd7ae7SLois Curfman McInnes {
88344cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense *) inA->data;
88444cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->A->data;
88544cd7ae7SLois Curfman McInnes   int          one = 1, nz;
88644cd7ae7SLois Curfman McInnes 
8873a40ed3dSBarry Smith   PetscFunctionBegin;
88844cd7ae7SLois Curfman McInnes   nz = a->m*a->n;
88944cd7ae7SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
89044cd7ae7SLois Curfman McInnes   PLogFlops(nz);
8913a40ed3dSBarry Smith   PetscFunctionReturn(0);
89244cd7ae7SLois Curfman McInnes }
89344cd7ae7SLois Curfman McInnes 
8943d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat,Mat *,int);
8958965ea79SLois Curfman McInnes 
8968965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
89709dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
89809dc0095SBarry Smith        MatGetRow_MPIDense,
89909dc0095SBarry Smith        MatRestoreRow_MPIDense,
90009dc0095SBarry Smith        MatMult_MPIDense,
90109dc0095SBarry Smith        MatMultAdd_MPIDense,
90209dc0095SBarry Smith        MatMultTrans_MPIDense,
90309dc0095SBarry Smith        MatMultTransAdd_MPIDense,
9048965ea79SLois Curfman McInnes        0,
90509dc0095SBarry Smith        0,
90609dc0095SBarry Smith        0,
90709dc0095SBarry Smith        0,
90809dc0095SBarry Smith        0,
90909dc0095SBarry Smith        0,
91009dc0095SBarry Smith        0,
91109dc0095SBarry Smith        MatTranspose_MPIDense,
91209dc0095SBarry Smith        MatGetInfo_MPIDense,0,
91309dc0095SBarry Smith        MatGetDiagonal_MPIDense,
91409dc0095SBarry Smith        0,
91509dc0095SBarry Smith        MatNorm_MPIDense,
91609dc0095SBarry Smith        MatAssemblyBegin_MPIDense,
91709dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
91809dc0095SBarry Smith        0,
91909dc0095SBarry Smith        MatSetOption_MPIDense,
92009dc0095SBarry Smith        MatZeroEntries_MPIDense,
92109dc0095SBarry Smith        MatZeroRows_MPIDense,
92209dc0095SBarry Smith        0,
92309dc0095SBarry Smith        0,
92409dc0095SBarry Smith        0,
92509dc0095SBarry Smith        0,
92609dc0095SBarry Smith        MatGetSize_MPIDense,
92709dc0095SBarry Smith        MatGetLocalSize_MPIDense,
92839ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
92909dc0095SBarry Smith        0,
93009dc0095SBarry Smith        0,
93109dc0095SBarry Smith        MatGetArray_MPIDense,
93209dc0095SBarry Smith        MatRestoreArray_MPIDense,
93394a9d846SBarry Smith        MatConvertSameType_MPIDense,
93409dc0095SBarry Smith        0,
93509dc0095SBarry Smith        0,
93609dc0095SBarry Smith        0,
93709dc0095SBarry Smith        0,
93809dc0095SBarry Smith        0,
93909dc0095SBarry Smith        0,
94009dc0095SBarry Smith        0,
94109dc0095SBarry Smith        MatGetValues_MPIDense,
94209dc0095SBarry Smith        0,
94309dc0095SBarry Smith        0,
94409dc0095SBarry Smith        MatScale_MPIDense,
94509dc0095SBarry Smith        0,
94609dc0095SBarry Smith        0,
94709dc0095SBarry Smith        0,
94809dc0095SBarry Smith        MatGetBlockSize_MPIDense,
94909dc0095SBarry Smith        0,
95009dc0095SBarry Smith        0,
95109dc0095SBarry Smith        0,
95209dc0095SBarry Smith        0,
95309dc0095SBarry Smith        0,
95409dc0095SBarry Smith        0,
95509dc0095SBarry Smith        0,
95609dc0095SBarry Smith        0,
95709dc0095SBarry Smith        0,
95809dc0095SBarry Smith        0,
95909dc0095SBarry Smith        0,
96009dc0095SBarry Smith        0,
96109dc0095SBarry Smith        MatGetMaps_Petsc};
9628965ea79SLois Curfman McInnes 
9635615d1e5SSatish Balay #undef __FUNC__
9645615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense"
9658965ea79SLois Curfman McInnes /*@C
96639ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
9678965ea79SLois Curfman McInnes 
968db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
969db81eaa0SLois Curfman McInnes 
9708965ea79SLois Curfman McInnes    Input Parameters:
971db81eaa0SLois Curfman McInnes +  comm - MPI communicator
9728965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
973db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
9748965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
975db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
976db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
977dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
9788965ea79SLois Curfman McInnes 
9798965ea79SLois Curfman McInnes    Output Parameter:
980477f1c0bSLois Curfman McInnes .  A - the matrix
9818965ea79SLois Curfman McInnes 
982b259b22eSLois Curfman McInnes    Notes:
98339ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
98439ddd567SLois Curfman McInnes    storage by columns.
9858965ea79SLois Curfman McInnes 
98618f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
98718f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
988b4fd4287SBarry Smith    set data=PETSC_NULL.
98918f449edSLois Curfman McInnes 
9908965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
9918965ea79SLois Curfman McInnes    (possibly both).
9928965ea79SLois Curfman McInnes 
9933501a2bdSLois Curfman McInnes    Currently, the only parallel dense matrix decomposition is by rows,
9943501a2bdSLois Curfman McInnes    so that n=N and each submatrix owns all of the global columns.
9953501a2bdSLois Curfman McInnes 
99639ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel
9978965ea79SLois Curfman McInnes 
99839ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
9998965ea79SLois Curfman McInnes @*/
1000477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
10018965ea79SLois Curfman McInnes {
10028965ea79SLois Curfman McInnes   Mat          mat;
100339ddd567SLois Curfman McInnes   Mat_MPIDense *a;
100425cdf11fSBarry Smith   int          ierr, i,flg;
10058965ea79SLois Curfman McInnes 
10063a40ed3dSBarry Smith   PetscFunctionBegin;
1007ed2daf61SLois Curfman McInnes   /* Note:  For now, when data is specified above, this assumes the user correctly
1008ed2daf61SLois Curfman McInnes    allocates the local dense storage space.  We should add error checking. */
100918f449edSLois Curfman McInnes 
1010477f1c0bSLois Curfman McInnes   *A = 0;
1011f830108cSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,comm,MatDestroy,MatView);
10128965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
10130452661fSBarry Smith   mat->data       = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
101409dc0095SBarry Smith   PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));
1015e1311b90SBarry Smith   mat->ops->destroy    = MatDestroy_MPIDense;
1016e1311b90SBarry Smith   mat->ops->view       = MatView_MPIDense;
10178965ea79SLois Curfman McInnes   mat->factor          = 0;
101890f02eecSBarry Smith   mat->mapping         = 0;
10198965ea79SLois Curfman McInnes 
1020622d7880SLois Curfman McInnes   a->factor       = 0;
1021e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
10228965ea79SLois Curfman McInnes   MPI_Comm_rank(comm,&a->rank);
10238965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&a->size);
10248965ea79SLois Curfman McInnes 
1025ca161407SBarry Smith   if (M == PETSC_DECIDE) {ierr = MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);}
10268965ea79SLois Curfman McInnes   if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);}
102739ddd567SLois Curfman McInnes 
102839ddd567SLois Curfman McInnes   /* each row stores all columns */
102939ddd567SLois Curfman McInnes   if (N == PETSC_DECIDE) N = n;
1030d7e8b826SBarry Smith   if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);}
1031a8c6a408SBarry Smith   /*  if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */
1032aca0ad90SLois Curfman McInnes   a->N = mat->N = N;
1033aca0ad90SLois Curfman McInnes   a->M = mat->M = M;
1034aca0ad90SLois Curfman McInnes   a->m = mat->m = m;
1035aca0ad90SLois Curfman McInnes   a->n = mat->n = n;
10368965ea79SLois Curfman McInnes 
10378965ea79SLois Curfman McInnes   /* build local table of row and column ownerships */
1038d7e8b826SBarry Smith   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners);
1039d7e8b826SBarry Smith   a->cowners = a->rowners + a->size + 1;
1040f09e8eb9SSatish Balay   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1041ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
10428965ea79SLois Curfman McInnes   a->rowners[0] = 0;
10438965ea79SLois Curfman McInnes   for ( i=2; i<=a->size; i++ ) {
10448965ea79SLois Curfman McInnes     a->rowners[i] += a->rowners[i-1];
10458965ea79SLois Curfman McInnes   }
10468965ea79SLois Curfman McInnes   a->rstart = a->rowners[a->rank];
10478965ea79SLois Curfman McInnes   a->rend   = a->rowners[a->rank+1];
1048ca161407SBarry Smith   ierr      = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
1049d7e8b826SBarry Smith   a->cowners[0] = 0;
1050d7e8b826SBarry Smith   for ( i=2; i<=a->size; i++ ) {
1051d7e8b826SBarry Smith     a->cowners[i] += a->cowners[i-1];
1052d7e8b826SBarry Smith   }
10538965ea79SLois Curfman McInnes 
1054029af93fSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A); CHKERRQ(ierr);
10558965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
10568965ea79SLois Curfman McInnes 
10578965ea79SLois Curfman McInnes   /* build cache for off array entries formed */
10588965ea79SLois Curfman McInnes   ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr);
10598965ea79SLois Curfman McInnes 
10608965ea79SLois Curfman McInnes   /* stuff used for matrix vector multiply */
10618965ea79SLois Curfman McInnes   a->lvec        = 0;
10628965ea79SLois Curfman McInnes   a->Mvctx       = 0;
106339b7565bSBarry Smith   a->roworiented = 1;
10648965ea79SLois Curfman McInnes 
1065477f1c0bSLois Curfman McInnes   *A = mat;
106625cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
106725cdf11fSBarry Smith   if (flg) {
10688c469469SLois Curfman McInnes     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
10698c469469SLois Curfman McInnes   }
10703a40ed3dSBarry Smith   PetscFunctionReturn(0);
10718965ea79SLois Curfman McInnes }
10728965ea79SLois Curfman McInnes 
10735615d1e5SSatish Balay #undef __FUNC__
10745615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIDense"
10753d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat A,Mat *newmat,int cpvalues)
10768965ea79SLois Curfman McInnes {
10778965ea79SLois Curfman McInnes   Mat          mat;
10783501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data;
107939ddd567SLois Curfman McInnes   int          ierr;
10802ba99913SLois Curfman McInnes   FactorCtx    *factor;
10818965ea79SLois Curfman McInnes 
10823a40ed3dSBarry Smith   PetscFunctionBegin;
10838965ea79SLois Curfman McInnes   *newmat       = 0;
1084f830108cSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,A->comm,MatDestroy,MatView);
10858965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
10860452661fSBarry Smith   mat->data      = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
108709dc0095SBarry Smith   PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));
1088e1311b90SBarry Smith   mat->ops->destroy   = MatDestroy_MPIDense;
1089e1311b90SBarry Smith   mat->ops->view      = MatView_MPIDense;
10903501a2bdSLois Curfman McInnes   mat->factor         = A->factor;
1091c456f294SBarry Smith   mat->assembled      = PETSC_TRUE;
10928965ea79SLois Curfman McInnes 
109344cd7ae7SLois Curfman McInnes   a->m = mat->m = oldmat->m;
109444cd7ae7SLois Curfman McInnes   a->n = mat->n = oldmat->n;
109544cd7ae7SLois Curfman McInnes   a->M = mat->M = oldmat->M;
109644cd7ae7SLois Curfman McInnes   a->N = mat->N = oldmat->N;
10972ba99913SLois Curfman McInnes   if (oldmat->factor) {
10982ba99913SLois Curfman McInnes     a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx)); CHKPTRQ(factor);
10992ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
11002ba99913SLois Curfman McInnes   } else a->factor = 0;
11018965ea79SLois Curfman McInnes 
11028965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
11038965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
11048965ea79SLois Curfman McInnes   a->size         = oldmat->size;
11058965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1106e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
11078965ea79SLois Curfman McInnes 
11080452661fSBarry Smith   a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners);
1109f09e8eb9SSatish Balay   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
11108965ea79SLois Curfman McInnes   PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));
11118965ea79SLois Curfman McInnes   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
11128965ea79SLois Curfman McInnes 
11138965ea79SLois Curfman McInnes   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
11148965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->lvec);
111555659b69SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
11168965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->Mvctx);
11178965ea79SLois Curfman McInnes   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
11188965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
11198965ea79SLois Curfman McInnes   *newmat = mat;
11203a40ed3dSBarry Smith   PetscFunctionReturn(0);
11218965ea79SLois Curfman McInnes }
11228965ea79SLois Curfman McInnes 
112377c4ece6SBarry Smith #include "sys.h"
11248965ea79SLois Curfman McInnes 
11255615d1e5SSatish Balay #undef __FUNC__
11265615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile"
112790ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat)
112890ace30eSBarry Smith {
112940011551SBarry Smith   int        *rowners, i,size,rank,m,ierr,nz,j;
113090ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
113190ace30eSBarry Smith   MPI_Status status;
113290ace30eSBarry Smith 
11333a40ed3dSBarry Smith   PetscFunctionBegin;
113490ace30eSBarry Smith   MPI_Comm_rank(comm,&rank);
113590ace30eSBarry Smith   MPI_Comm_size(comm,&size);
113690ace30eSBarry Smith 
113790ace30eSBarry Smith   /* determine ownership of all rows */
113890ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
113990ace30eSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
1140ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
114190ace30eSBarry Smith   rowners[0] = 0;
114290ace30eSBarry Smith   for ( i=2; i<=size; i++ ) {
114390ace30eSBarry Smith     rowners[i] += rowners[i-1];
114490ace30eSBarry Smith   }
114590ace30eSBarry Smith 
114690ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
114790ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array); CHKERRQ(ierr);
114890ace30eSBarry Smith 
114990ace30eSBarry Smith   if (!rank) {
115090ace30eSBarry Smith     vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
115190ace30eSBarry Smith 
115290ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
11530752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR); CHKERRQ(ierr);
115490ace30eSBarry Smith 
115590ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
115690ace30eSBarry Smith     vals_ptr = vals;
115790ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
115890ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
115990ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
116090ace30eSBarry Smith       }
116190ace30eSBarry Smith     }
116290ace30eSBarry Smith 
116390ace30eSBarry Smith     /* read in other processors and ship out */
116490ace30eSBarry Smith     for ( i=1; i<size; i++ ) {
116590ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
11660752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
1167ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
116890ace30eSBarry Smith     }
11693a40ed3dSBarry Smith   } else {
117090ace30eSBarry Smith     /* receive numeric values */
117190ace30eSBarry Smith     vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
117290ace30eSBarry Smith 
117390ace30eSBarry Smith     /* receive message of values*/
1174ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
117590ace30eSBarry Smith 
117690ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
117790ace30eSBarry Smith     vals_ptr = vals;
117890ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
117990ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
118090ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
118190ace30eSBarry Smith       }
118290ace30eSBarry Smith     }
118390ace30eSBarry Smith   }
118490ace30eSBarry Smith   PetscFree(rowners);
118590ace30eSBarry Smith   PetscFree(vals);
11866d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11876d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11883a40ed3dSBarry Smith   PetscFunctionReturn(0);
118990ace30eSBarry Smith }
119090ace30eSBarry Smith 
119190ace30eSBarry Smith 
11925615d1e5SSatish Balay #undef __FUNC__
11935615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense"
119419bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
11958965ea79SLois Curfman McInnes {
11968965ea79SLois Curfman McInnes   Mat          A;
11978965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
119819bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
11998965ea79SLois Curfman McInnes   MPI_Status   status;
12008965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
12018965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
120219bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
12033a40ed3dSBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
12048965ea79SLois Curfman McInnes 
12053a40ed3dSBarry Smith   PetscFunctionBegin;
12068965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank);
12078965ea79SLois Curfman McInnes   if (!rank) {
120890ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
12090752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr);
1210a8c6a408SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object");
12118965ea79SLois Curfman McInnes   }
12128965ea79SLois Curfman McInnes 
1213ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
121490ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
121590ace30eSBarry Smith 
121690ace30eSBarry Smith   /*
121790ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
121890ace30eSBarry Smith   */
121990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
12203a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
12213a40ed3dSBarry Smith     PetscFunctionReturn(0);
122290ace30eSBarry Smith   }
122390ace30eSBarry Smith 
12248965ea79SLois Curfman McInnes   /* determine ownership of all rows */
12258965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
12260452661fSBarry Smith   rowners    = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
1227ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
12288965ea79SLois Curfman McInnes   rowners[0] = 0;
12298965ea79SLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
12308965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
12318965ea79SLois Curfman McInnes   }
12328965ea79SLois Curfman McInnes   rstart = rowners[rank];
12338965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
12348965ea79SLois Curfman McInnes 
12358965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
12360452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
12378965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
12388965ea79SLois Curfman McInnes   if (!rank) {
12390452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
12400752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
12410452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts);
12428965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1243ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
12440452661fSBarry Smith     PetscFree(sndcounts);
1245ca161407SBarry Smith   } else {
1246ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr);
12478965ea79SLois Curfman McInnes   }
12488965ea79SLois Curfman McInnes 
12498965ea79SLois Curfman McInnes   if (!rank) {
12508965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
12510452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz);
1252cddf8d76SBarry Smith     PetscMemzero(procsnz,size*sizeof(int));
12538965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
12548965ea79SLois Curfman McInnes       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
12558965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
12568965ea79SLois Curfman McInnes       }
12578965ea79SLois Curfman McInnes     }
12580452661fSBarry Smith     PetscFree(rowlengths);
12598965ea79SLois Curfman McInnes 
12608965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
12618965ea79SLois Curfman McInnes     maxnz = 0;
12628965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
12630452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
12648965ea79SLois Curfman McInnes     }
12650452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols);
12668965ea79SLois Curfman McInnes 
12678965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
12688965ea79SLois Curfman McInnes     nz = procsnz[0];
12690452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
12700752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr);
12718965ea79SLois Curfman McInnes 
12728965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
12738965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
12748965ea79SLois Curfman McInnes       nz   = procsnz[i];
12750752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
1276ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
12778965ea79SLois Curfman McInnes     }
12780452661fSBarry Smith     PetscFree(cols);
12793a40ed3dSBarry Smith   } else {
12808965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
12818965ea79SLois Curfman McInnes     nz = 0;
12828965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
12838965ea79SLois Curfman McInnes       nz += ourlens[i];
12848965ea79SLois Curfman McInnes     }
12850452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
12868965ea79SLois Curfman McInnes 
12878965ea79SLois Curfman McInnes     /* receive message of column indices*/
1288ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1289ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
1290a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
12918965ea79SLois Curfman McInnes   }
12928965ea79SLois Curfman McInnes 
12938965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1294cddf8d76SBarry Smith   PetscMemzero(offlens,m*sizeof(int));
12958965ea79SLois Curfman McInnes   jj = 0;
12968965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
12978965ea79SLois Curfman McInnes     for ( j=0; j<ourlens[i]; j++ ) {
12988965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
12998965ea79SLois Curfman McInnes       jj++;
13008965ea79SLois Curfman McInnes     }
13018965ea79SLois Curfman McInnes   }
13028965ea79SLois Curfman McInnes 
13038965ea79SLois Curfman McInnes   /* create our matrix */
13048965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13058965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
13068965ea79SLois Curfman McInnes   }
1307b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
13088965ea79SLois Curfman McInnes   A = *newmat;
13098965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
13108965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
13118965ea79SLois Curfman McInnes   }
13128965ea79SLois Curfman McInnes 
13138965ea79SLois Curfman McInnes   if (!rank) {
13140452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
13158965ea79SLois Curfman McInnes 
13168965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
13178965ea79SLois Curfman McInnes     nz = procsnz[0];
13180752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
13198965ea79SLois Curfman McInnes 
13208965ea79SLois Curfman McInnes     /* insert into matrix */
13218965ea79SLois Curfman McInnes     jj      = rstart;
13228965ea79SLois Curfman McInnes     smycols = mycols;
13238965ea79SLois Curfman McInnes     svals   = vals;
13248965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
13258965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13268965ea79SLois Curfman McInnes       smycols += ourlens[i];
13278965ea79SLois Curfman McInnes       svals   += ourlens[i];
13288965ea79SLois Curfman McInnes       jj++;
13298965ea79SLois Curfman McInnes     }
13308965ea79SLois Curfman McInnes 
13318965ea79SLois Curfman McInnes     /* read in other processors and ship out */
13328965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
13338965ea79SLois Curfman McInnes       nz   = procsnz[i];
13340752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
1335ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
13368965ea79SLois Curfman McInnes     }
13370452661fSBarry Smith     PetscFree(procsnz);
13383a40ed3dSBarry Smith   } else {
13398965ea79SLois Curfman McInnes     /* receive numeric values */
13400452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
13418965ea79SLois Curfman McInnes 
13428965ea79SLois Curfman McInnes     /* receive message of values*/
1343ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1344ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
1345a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
13468965ea79SLois Curfman McInnes 
13478965ea79SLois Curfman McInnes     /* insert into matrix */
13488965ea79SLois Curfman McInnes     jj      = rstart;
13498965ea79SLois Curfman McInnes     smycols = mycols;
13508965ea79SLois Curfman McInnes     svals   = vals;
13518965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
13528965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
13538965ea79SLois Curfman McInnes       smycols += ourlens[i];
13548965ea79SLois Curfman McInnes       svals   += ourlens[i];
13558965ea79SLois Curfman McInnes       jj++;
13568965ea79SLois Curfman McInnes     }
13578965ea79SLois Curfman McInnes   }
13580452661fSBarry Smith   PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners);
13598965ea79SLois Curfman McInnes 
13606d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13616d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13623a40ed3dSBarry Smith   PetscFunctionReturn(0);
13638965ea79SLois Curfman McInnes }
136490ace30eSBarry Smith 
136590ace30eSBarry Smith 
136690ace30eSBarry Smith 
136790ace30eSBarry Smith 
136890ace30eSBarry Smith 
1369