xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 639f9d9dbbc54d6ac4e42e98283c540b41bb2cee)
18965ea79SLois Curfman McInnes #ifndef lint
2*639f9d9dSBarry Smith static char vcid[] = "$Id: mpidense.c,v 1.50 1996/10/01 03:34:29 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 
127056b6fcSBarry Smith static int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
138965ea79SLois Curfman McInnes {
1439b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense *) mat->data;
1539b7565bSBarry Smith   int          ierr, i, j, rstart = A->rstart, rend = A->rend, row;
1639b7565bSBarry Smith   int          roworiented = A->roworiented;
178965ea79SLois Curfman McInnes 
1839b7565bSBarry Smith   if (A->insertmode != NOT_SET_VALUES && A->insertmode != addv) {
1939ddd567SLois Curfman McInnes     SETERRQ(1,"MatSetValues_MPIDense:Cannot mix inserts and adds");
208965ea79SLois Curfman McInnes   }
2139b7565bSBarry Smith   A->insertmode = addv;
228965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
2339ddd567SLois Curfman McInnes     if (idxm[i] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative row");
2439b7565bSBarry Smith     if (idxm[i] >= A->M) SETERRQ(1,"MatSetValues_MPIDense:Row too large");
258965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
268965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
2739b7565bSBarry Smith       if (roworiented) {
2839b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv); CHKERRQ(ierr);
2939b7565bSBarry Smith       }
3039b7565bSBarry Smith       else {
318965ea79SLois Curfman McInnes         for ( j=0; j<n; j++ ) {
3239ddd567SLois Curfman McInnes           if (idxn[j] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative column");
3339b7565bSBarry Smith           if (idxn[j] >= A->N) SETERRQ(1,"MatSetValues_MPIDense:Column too large");
3439b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv); CHKERRQ(ierr);
3539b7565bSBarry Smith         }
368965ea79SLois Curfman McInnes       }
378965ea79SLois Curfman McInnes     }
388965ea79SLois Curfman McInnes     else {
3939b7565bSBarry Smith       if (roworiented) {
4039b7565bSBarry Smith         ierr = StashValues_Private(&A->stash,idxm[i],n,idxn,v+i*n,addv); CHKERRQ(ierr);
4139b7565bSBarry Smith       }
4239b7565bSBarry Smith       else { /* must stash each seperately */
4339b7565bSBarry Smith         row = idxm[i];
4439b7565bSBarry Smith         for ( j=0; j<n; j++ ) {
457056b6fcSBarry Smith           ierr = StashValues_Private(&A->stash,row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
4639b7565bSBarry Smith         }
4739b7565bSBarry Smith       }
48b49de8d1SLois Curfman McInnes     }
49b49de8d1SLois Curfman McInnes   }
50b49de8d1SLois Curfman McInnes   return 0;
51b49de8d1SLois Curfman McInnes }
52b49de8d1SLois Curfman McInnes 
53b49de8d1SLois Curfman McInnes static int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
54b49de8d1SLois Curfman McInnes {
55b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
56b49de8d1SLois Curfman McInnes   int          ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row;
57b49de8d1SLois Curfman McInnes 
58b49de8d1SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
59b49de8d1SLois Curfman McInnes     if (idxm[i] < 0) SETERRQ(1,"MatGetValues_MPIDense:Negative row");
60b49de8d1SLois Curfman McInnes     if (idxm[i] >= mdn->M) SETERRQ(1,"MatGetValues_MPIDense:Row too large");
61b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
62b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
63b49de8d1SLois Curfman McInnes       for ( j=0; j<n; j++ ) {
64b49de8d1SLois Curfman McInnes         if (idxn[j] < 0) SETERRQ(1,"MatGetValues_MPIDense:Negative column");
65b49de8d1SLois Curfman McInnes         if (idxn[j] >= mdn->N)
66b49de8d1SLois Curfman McInnes           SETERRQ(1,"MatGetValues_MPIDense:Column too large");
67b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j); CHKERRQ(ierr);
68b49de8d1SLois Curfman McInnes       }
69b49de8d1SLois Curfman McInnes     }
70b49de8d1SLois Curfman McInnes     else {
71b49de8d1SLois Curfman McInnes       SETERRQ(1,"MatGetValues_MPIDense:Only local values currently supported");
728965ea79SLois Curfman McInnes     }
738965ea79SLois Curfman McInnes   }
748965ea79SLois Curfman McInnes   return 0;
758965ea79SLois Curfman McInnes }
768965ea79SLois Curfman McInnes 
77ff14e315SSatish Balay static int MatGetArray_MPIDense(Mat A,Scalar **array)
78ff14e315SSatish Balay {
79ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
80ff14e315SSatish Balay   int ierr;
81ff14e315SSatish Balay 
82ff14e315SSatish Balay   ierr = MatGetArray(a->A,array); CHKERRQ(ierr);
83ff14e315SSatish Balay   return 0;
84ff14e315SSatish Balay }
85ff14e315SSatish Balay 
86ff14e315SSatish Balay static int MatRestoreArray_MPIDense(Mat A,Scalar **array)
87ff14e315SSatish Balay {
88ff14e315SSatish Balay   return 0;
89ff14e315SSatish Balay }
90ff14e315SSatish Balay 
9139ddd567SLois Curfman McInnes static int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
928965ea79SLois Curfman McInnes {
9339ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
948965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
9539ddd567SLois Curfman McInnes   int          size = mdn->size, *owners = mdn->rowners, rank = mdn->rank;
968965ea79SLois Curfman McInnes   int          *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
9739ddd567SLois Curfman McInnes   int          tag = mat->tag, *owner,*starts,count,ierr;
988965ea79SLois Curfman McInnes   InsertMode   addv;
9939ddd567SLois Curfman McInnes   MPI_Request  *send_waits,*recv_waits;
1008965ea79SLois Curfman McInnes   Scalar       *rvalues,*svalues;
1018965ea79SLois Curfman McInnes 
1028965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
1036d84be18SBarry Smith   MPI_Allreduce(&mdn->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
1047056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
1057056b6fcSBarry Smith     SETERRQ(1,"MatAssemblyBegin_MPIDense:Cannot mix adds/inserts on different procs");
1068965ea79SLois Curfman McInnes   }
10739ddd567SLois Curfman McInnes   mdn->insertmode = addv; /* in case this processor had no cache */
1088965ea79SLois Curfman McInnes 
1098965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
1100452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
111cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
1120452661fSBarry Smith   owner = (int *) PetscMalloc( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
11339ddd567SLois Curfman McInnes   for ( i=0; i<mdn->stash.n; i++ ) {
11439ddd567SLois Curfman McInnes     idx = mdn->stash.idx[i];
1158965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
1168965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
1178965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
1188965ea79SLois Curfman McInnes       }
1198965ea79SLois Curfman McInnes     }
1208965ea79SLois Curfman McInnes   }
1218965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
1228965ea79SLois Curfman McInnes 
1238965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
1240452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
1256d84be18SBarry Smith   MPI_Allreduce(procs,work,size,MPI_INT,MPI_SUM,comm);
1268965ea79SLois Curfman McInnes   nreceives = work[rank];
1273b2fbd54SBarry Smith   if (nreceives > size) SETERRQ(1,"MatAssemblyBegin_MPIDense:Internal PETSc error");
1286d84be18SBarry Smith   MPI_Allreduce(nprocs,work,size,MPI_INT,MPI_MAX,comm);
1298965ea79SLois Curfman McInnes   nmax = work[rank];
1300452661fSBarry Smith   PetscFree(work);
1318965ea79SLois Curfman McInnes 
1328965ea79SLois Curfman McInnes   /* post receives:
1338965ea79SLois Curfman McInnes        1) each message will consist of ordered pairs
1348965ea79SLois Curfman McInnes      (global index,value) we store the global index as a double
1358965ea79SLois Curfman McInnes      to simplify the message passing.
1368965ea79SLois Curfman McInnes        2) since we don't know how long each individual message is we
1378965ea79SLois Curfman McInnes      allocate the largest needed buffer for each receive. Potentially
1388965ea79SLois Curfman McInnes      this is a lot of wasted space.
1398965ea79SLois Curfman McInnes 
1408965ea79SLois Curfman McInnes        This could be done better.
1418965ea79SLois Curfman McInnes   */
1423b2fbd54SBarry Smith   rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));CHKPTRQ(rvalues);
1433b2fbd54SBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
1448965ea79SLois Curfman McInnes   for ( i=0; i<nreceives; i++ ) {
1457056b6fcSBarry Smith     MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
1468965ea79SLois Curfman McInnes   }
1478965ea79SLois Curfman McInnes 
1488965ea79SLois Curfman McInnes   /* do sends:
1498965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
1508965ea79SLois Curfman McInnes          the ith processor
1518965ea79SLois Curfman McInnes   */
1523b2fbd54SBarry Smith   svalues = (Scalar *) PetscMalloc( 3*(mdn->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues);
1533b2fbd54SBarry Smith   send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
1540452661fSBarry Smith   starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts);
1558965ea79SLois Curfman McInnes   starts[0] = 0;
1568965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
15739ddd567SLois Curfman McInnes   for ( i=0; i<mdn->stash.n; i++ ) {
15839ddd567SLois Curfman McInnes     svalues[3*starts[owner[i]]]       = (Scalar)  mdn->stash.idx[i];
15939ddd567SLois Curfman McInnes     svalues[3*starts[owner[i]]+1]     = (Scalar)  mdn->stash.idy[i];
16039ddd567SLois Curfman McInnes     svalues[3*(starts[owner[i]]++)+2] =  mdn->stash.array[i];
1618965ea79SLois Curfman McInnes   }
1620452661fSBarry Smith   PetscFree(owner);
1638965ea79SLois Curfman McInnes   starts[0] = 0;
1648965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
1658965ea79SLois Curfman McInnes   count = 0;
1668965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
1678965ea79SLois Curfman McInnes     if (procs[i]) {
1687056b6fcSBarry Smith       MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);
1698965ea79SLois Curfman McInnes     }
1708965ea79SLois Curfman McInnes   }
1710452661fSBarry Smith   PetscFree(starts); PetscFree(nprocs);
1728965ea79SLois Curfman McInnes 
1738965ea79SLois Curfman McInnes   /* Free cache space */
17439ddd567SLois Curfman McInnes   ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr);
1758965ea79SLois Curfman McInnes 
17639ddd567SLois Curfman McInnes   mdn->svalues    = svalues;    mdn->rvalues = rvalues;
17739ddd567SLois Curfman McInnes   mdn->nsends     = nsends;     mdn->nrecvs = nreceives;
17839ddd567SLois Curfman McInnes   mdn->send_waits = send_waits; mdn->recv_waits = recv_waits;
17939ddd567SLois Curfman McInnes   mdn->rmax       = nmax;
1808965ea79SLois Curfman McInnes 
1818965ea79SLois Curfman McInnes   return 0;
1828965ea79SLois Curfman McInnes }
18339ddd567SLois Curfman McInnes extern int MatSetUpMultiply_MPIDense(Mat);
1848965ea79SLois Curfman McInnes 
18539ddd567SLois Curfman McInnes static int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
1868965ea79SLois Curfman McInnes {
18739ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
1888965ea79SLois Curfman McInnes   MPI_Status   *send_status,recv_status;
18939ddd567SLois Curfman McInnes   int          imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col;
1908965ea79SLois Curfman McInnes   Scalar       *values,val;
19139ddd567SLois Curfman McInnes   InsertMode   addv = mdn->insertmode;
1928965ea79SLois Curfman McInnes 
1938965ea79SLois Curfman McInnes   /*  wait on receives */
1948965ea79SLois Curfman McInnes   while (count) {
19539ddd567SLois Curfman McInnes     MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status);
1968965ea79SLois Curfman McInnes     /* unpack receives into our local space */
19739ddd567SLois Curfman McInnes     values = mdn->rvalues + 3*imdex*mdn->rmax;
1988965ea79SLois Curfman McInnes     MPI_Get_count(&recv_status,MPIU_SCALAR,&n);
1998965ea79SLois Curfman McInnes     n = n/3;
2008965ea79SLois Curfman McInnes     for ( i=0; i<n; i++ ) {
201227d817aSBarry Smith       row = (int) PetscReal(values[3*i]) - mdn->rstart;
202227d817aSBarry Smith       col = (int) PetscReal(values[3*i+1]);
2038965ea79SLois Curfman McInnes       val = values[3*i+2];
20439ddd567SLois Curfman McInnes       if (col >= 0 && col < mdn->N) {
20539ddd567SLois Curfman McInnes         MatSetValues(mdn->A,1,&row,1,&col,&val,addv);
2068965ea79SLois Curfman McInnes       }
20739ddd567SLois Curfman McInnes       else {SETERRQ(1,"MatAssemblyEnd_MPIDense:Invalid column");}
2088965ea79SLois Curfman McInnes     }
2098965ea79SLois Curfman McInnes     count--;
2108965ea79SLois Curfman McInnes   }
2110452661fSBarry Smith   PetscFree(mdn->recv_waits); PetscFree(mdn->rvalues);
2128965ea79SLois Curfman McInnes 
2138965ea79SLois Curfman McInnes   /* wait on sends */
21439ddd567SLois Curfman McInnes   if (mdn->nsends) {
2157056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(mdn->nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
21639ddd567SLois Curfman McInnes     MPI_Waitall(mdn->nsends,mdn->send_waits,send_status);
2170452661fSBarry Smith     PetscFree(send_status);
2188965ea79SLois Curfman McInnes   }
2190452661fSBarry Smith   PetscFree(mdn->send_waits); PetscFree(mdn->svalues);
2208965ea79SLois Curfman McInnes 
22139ddd567SLois Curfman McInnes   mdn->insertmode = NOT_SET_VALUES;
22239ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr);
22339ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr);
2248965ea79SLois Curfman McInnes 
2256d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
22639ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr);
2278965ea79SLois Curfman McInnes   }
2288965ea79SLois Curfman McInnes   return 0;
2298965ea79SLois Curfman McInnes }
2308965ea79SLois Curfman McInnes 
23139ddd567SLois Curfman McInnes static int MatZeroEntries_MPIDense(Mat A)
2328965ea79SLois Curfman McInnes {
23339ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense *) A->data;
23439ddd567SLois Curfman McInnes   return MatZeroEntries(l->A);
2358965ea79SLois Curfman McInnes }
2368965ea79SLois Curfman McInnes 
2374e220ebcSLois Curfman McInnes static int MatGetBlockSize_MPIDense(Mat A,int *bs)
2384e220ebcSLois Curfman McInnes {
2394e220ebcSLois Curfman McInnes   *bs = 1;
2404e220ebcSLois Curfman McInnes   return 0;
2414e220ebcSLois Curfman McInnes }
2424e220ebcSLois Curfman McInnes 
2438965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2448965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2458965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2463501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2478965ea79SLois Curfman McInnes    routine.
2488965ea79SLois Curfman McInnes */
24939ddd567SLois Curfman McInnes static int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag)
2508965ea79SLois Curfman McInnes {
25139ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense *) A->data;
2528965ea79SLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
2538965ea79SLois Curfman McInnes   int            *procs,*nprocs,j,found,idx,nsends,*work;
2548965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2558965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2568965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2578965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2588965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2598965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2608965ea79SLois Curfman McInnes   IS             istmp;
2618965ea79SLois Curfman McInnes 
26277c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
2638965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
2648965ea79SLois Curfman McInnes 
2658965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
2660452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
267cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
2680452661fSBarry Smith   owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
2698965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
2708965ea79SLois Curfman McInnes     idx = rows[i];
2718965ea79SLois Curfman McInnes     found = 0;
2728965ea79SLois Curfman McInnes     for ( j=0; j<size; j++ ) {
2738965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
2748965ea79SLois Curfman McInnes         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
2758965ea79SLois Curfman McInnes       }
2768965ea79SLois Curfman McInnes     }
27739ddd567SLois Curfman McInnes     if (!found) SETERRQ(1,"MatZeroRows_MPIDense:Index out of range");
2788965ea79SLois Curfman McInnes   }
2798965ea79SLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
2808965ea79SLois Curfman McInnes 
2818965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
2820452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
2838965ea79SLois Curfman McInnes   MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);
2848965ea79SLois Curfman McInnes   nrecvs = work[rank];
2858965ea79SLois Curfman McInnes   MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);
2868965ea79SLois Curfman McInnes   nmax = work[rank];
2870452661fSBarry Smith   PetscFree(work);
2888965ea79SLois Curfman McInnes 
2898965ea79SLois Curfman McInnes   /* post receives:   */
2900452661fSBarry Smith   rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */
2918965ea79SLois Curfman McInnes   CHKPTRQ(rvalues);
2920452661fSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));
2938965ea79SLois Curfman McInnes   CHKPTRQ(recv_waits);
2948965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
2958965ea79SLois Curfman McInnes     MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
2968965ea79SLois Curfman McInnes   }
2978965ea79SLois Curfman McInnes 
2988965ea79SLois Curfman McInnes   /* do sends:
2998965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3008965ea79SLois Curfman McInnes          the ith processor
3018965ea79SLois Curfman McInnes   */
3020452661fSBarry Smith   svalues    = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
3037056b6fcSBarry Smith   send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
3040452661fSBarry Smith   starts     = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts);
3058965ea79SLois Curfman McInnes   starts[0]  = 0;
3068965ea79SLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3078965ea79SLois Curfman McInnes   for ( i=0; i<N; i++ ) {
3088965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3098965ea79SLois Curfman McInnes   }
3108965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3118965ea79SLois Curfman McInnes 
3128965ea79SLois Curfman McInnes   starts[0] = 0;
3138965ea79SLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3148965ea79SLois Curfman McInnes   count = 0;
3158965ea79SLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3168965ea79SLois Curfman McInnes     if (procs[i]) {
3178965ea79SLois Curfman McInnes       MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);
3188965ea79SLois Curfman McInnes     }
3198965ea79SLois Curfman McInnes   }
3200452661fSBarry Smith   PetscFree(starts);
3218965ea79SLois Curfman McInnes 
3228965ea79SLois Curfman McInnes   base = owners[rank];
3238965ea79SLois Curfman McInnes 
3248965ea79SLois Curfman McInnes   /*  wait on receives */
3250452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
3268965ea79SLois Curfman McInnes   source = lens + nrecvs;
3278965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3288965ea79SLois Curfman McInnes   while (count) {
3298965ea79SLois Curfman McInnes     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
3308965ea79SLois Curfman McInnes     /* unpack receives into our local space */
3318965ea79SLois Curfman McInnes     MPI_Get_count(&recv_status,MPI_INT,&n);
3328965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3338965ea79SLois Curfman McInnes     lens[imdex]  = n;
3348965ea79SLois Curfman McInnes     slen += n;
3358965ea79SLois Curfman McInnes     count--;
3368965ea79SLois Curfman McInnes   }
3370452661fSBarry Smith   PetscFree(recv_waits);
3388965ea79SLois Curfman McInnes 
3398965ea79SLois Curfman McInnes   /* move the data into the send scatter */
3400452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
3418965ea79SLois Curfman McInnes   count = 0;
3428965ea79SLois Curfman McInnes   for ( i=0; i<nrecvs; i++ ) {
3438965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3448965ea79SLois Curfman McInnes     for ( j=0; j<lens[i]; j++ ) {
3458965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3468965ea79SLois Curfman McInnes     }
3478965ea79SLois Curfman McInnes   }
3480452661fSBarry Smith   PetscFree(rvalues); PetscFree(lens);
3490452661fSBarry Smith   PetscFree(owner); PetscFree(nprocs);
3508965ea79SLois Curfman McInnes 
3518965ea79SLois Curfman McInnes   /* actually zap the local rows */
352537820f0SBarry Smith   ierr = ISCreateGeneral(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
3538965ea79SLois Curfman McInnes   PLogObjectParent(A,istmp);
3540452661fSBarry Smith   PetscFree(lrows);
3558965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
3568965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp); CHKERRQ(ierr);
3578965ea79SLois Curfman McInnes 
3588965ea79SLois Curfman McInnes   /* wait on sends */
3598965ea79SLois Curfman McInnes   if (nsends) {
3607056b6fcSBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
3618965ea79SLois Curfman McInnes     MPI_Waitall(nsends,send_waits,send_status);
3620452661fSBarry Smith     PetscFree(send_status);
3638965ea79SLois Curfman McInnes   }
3640452661fSBarry Smith   PetscFree(send_waits); PetscFree(svalues);
3658965ea79SLois Curfman McInnes 
3668965ea79SLois Curfman McInnes   return 0;
3678965ea79SLois Curfman McInnes }
3688965ea79SLois Curfman McInnes 
36939ddd567SLois Curfman McInnes static int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3708965ea79SLois Curfman McInnes {
37139ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
3728965ea79SLois Curfman McInnes   int          ierr;
373c456f294SBarry Smith 
3747056b6fcSBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
3757056b6fcSBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
37644cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy); CHKERRQ(ierr);
3778965ea79SLois Curfman McInnes   return 0;
3788965ea79SLois Curfman McInnes }
3798965ea79SLois Curfman McInnes 
38039ddd567SLois Curfman McInnes static int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
3818965ea79SLois Curfman McInnes {
38239ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
3838965ea79SLois Curfman McInnes   int          ierr;
384c456f294SBarry Smith 
385c456f294SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
386c456f294SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr);
38744cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr);
3888965ea79SLois Curfman McInnes   return 0;
3898965ea79SLois Curfman McInnes }
3908965ea79SLois Curfman McInnes 
391096963f5SLois Curfman McInnes static int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy)
392096963f5SLois Curfman McInnes {
393096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
394096963f5SLois Curfman McInnes   int          ierr;
3953501a2bdSLois Curfman McInnes   Scalar       zero = 0.0;
396096963f5SLois Curfman McInnes 
3973501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy); CHKERRQ(ierr);
39844cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr);
399537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
400537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
401096963f5SLois Curfman McInnes   return 0;
402096963f5SLois Curfman McInnes }
403096963f5SLois Curfman McInnes 
404096963f5SLois Curfman McInnes static int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
405096963f5SLois Curfman McInnes {
406096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
407096963f5SLois Curfman McInnes   int          ierr;
408096963f5SLois Curfman McInnes 
4093501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz); CHKERRQ(ierr);
41044cd7ae7SLois Curfman McInnes   ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr);
411537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
412537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
413096963f5SLois Curfman McInnes   return 0;
414096963f5SLois Curfman McInnes }
415096963f5SLois Curfman McInnes 
41639ddd567SLois Curfman McInnes static int MatGetDiagonal_MPIDense(Mat A,Vec v)
4178965ea79SLois Curfman McInnes {
41839ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
419096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data;
42044cd7ae7SLois Curfman McInnes   int          ierr, len, i, n, m = a->m, radd;
42144cd7ae7SLois Curfman McInnes   Scalar       *x, zero = 0.0;
422ed3cc1f0SBarry Smith 
42344cd7ae7SLois Curfman McInnes   VecSet(&zero,v);
424096963f5SLois Curfman McInnes   ierr = VecGetArray(v,&x); CHKERRQ(ierr);
425096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n); CHKERRQ(ierr);
426096963f5SLois Curfman McInnes   if (n != a->M) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
42744cd7ae7SLois Curfman McInnes   len = PetscMin(aloc->m,aloc->n);
4287ddc982cSLois Curfman McInnes   radd = a->rstart*m;
42944cd7ae7SLois Curfman McInnes   for ( i=0; i<len; i++ ) {
430096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
431096963f5SLois Curfman McInnes   }
432096963f5SLois Curfman McInnes   return 0;
4338965ea79SLois Curfman McInnes }
4348965ea79SLois Curfman McInnes 
43539ddd567SLois Curfman McInnes static int MatDestroy_MPIDense(PetscObject obj)
4368965ea79SLois Curfman McInnes {
4378965ea79SLois Curfman McInnes   Mat          mat = (Mat) obj;
4383501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4398965ea79SLois Curfman McInnes   int          ierr;
440ed3cc1f0SBarry Smith 
4418965ea79SLois Curfman McInnes #if defined(PETSC_LOG)
4423501a2bdSLois Curfman McInnes   PLogObjectState(obj,"Rows=%d, Cols=%d",mdn->M,mdn->N);
4438965ea79SLois Curfman McInnes #endif
4440452661fSBarry Smith   PetscFree(mdn->rowners);
4453501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A); CHKERRQ(ierr);
4463501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4473501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
448622d7880SLois Curfman McInnes   if (mdn->factor) {
449622d7880SLois Curfman McInnes     if (mdn->factor->temp)   PetscFree(mdn->factor->temp);
450622d7880SLois Curfman McInnes     if (mdn->factor->tag)    PetscFree(mdn->factor->tag);
451622d7880SLois Curfman McInnes     if (mdn->factor->pivots) PetscFree(mdn->factor->pivots);
452622d7880SLois Curfman McInnes     PetscFree(mdn->factor);
453622d7880SLois Curfman McInnes   }
4540452661fSBarry Smith   PetscFree(mdn);
4558965ea79SLois Curfman McInnes   PLogObjectDestroy(mat);
4560452661fSBarry Smith   PetscHeaderDestroy(mat);
4578965ea79SLois Curfman McInnes   return 0;
4588965ea79SLois Curfman McInnes }
45939ddd567SLois Curfman McInnes 
4608965ea79SLois Curfman McInnes #include "pinclude/pviewer.h"
4618965ea79SLois Curfman McInnes 
46239ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer)
4638965ea79SLois Curfman McInnes {
46439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
4658965ea79SLois Curfman McInnes   int          ierr;
4667056b6fcSBarry Smith 
46739ddd567SLois Curfman McInnes   if (mdn->size == 1) {
46839ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
4698965ea79SLois Curfman McInnes   }
47039ddd567SLois Curfman McInnes   else SETERRQ(1,"MatView_MPIDense_Binary:Only uniprocessor output supported");
4718965ea79SLois Curfman McInnes   return 0;
4728965ea79SLois Curfman McInnes }
4738965ea79SLois Curfman McInnes 
47439ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer)
4758965ea79SLois Curfman McInnes {
47639ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
47739ddd567SLois Curfman McInnes   int          ierr, format;
4788965ea79SLois Curfman McInnes   FILE         *fd;
47919bcc07fSBarry Smith   ViewerType   vtype;
4808965ea79SLois Curfman McInnes 
48119bcc07fSBarry Smith   ViewerGetType(viewer,&vtype);
48290ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
48390ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
484*639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
4854e220ebcSLois Curfman McInnes     int rank;
4864e220ebcSLois Curfman McInnes     MatInfo info;
487096963f5SLois Curfman McInnes     MPI_Comm_rank(mat->comm,&rank);
4884e220ebcSLois Curfman McInnes     ierr = MatGetInfo(mat,MAT_LOCAL,&info);
48977c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
4904e220ebcSLois Curfman McInnes       fprintf(fd,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m,
4914e220ebcSLois Curfman McInnes          (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
492096963f5SLois Curfman McInnes       fflush(fd);
49377c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
4943501a2bdSLois Curfman McInnes     ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr);
4953501a2bdSLois Curfman McInnes     return 0;
4963501a2bdSLois Curfman McInnes   }
49702cad45dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_INFO) {
498096963f5SLois Curfman McInnes     return 0;
4998965ea79SLois Curfman McInnes   }
50019bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
50177c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
50239ddd567SLois Curfman McInnes     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n",
50339ddd567SLois Curfman McInnes              mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n);
50439ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5058965ea79SLois Curfman McInnes     fflush(fd);
50677c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
5078965ea79SLois Curfman McInnes   }
5088965ea79SLois Curfman McInnes   else {
50939ddd567SLois Curfman McInnes     int size = mdn->size, rank = mdn->rank;
5108965ea79SLois Curfman McInnes     if (size == 1) {
51139ddd567SLois Curfman McInnes       ierr = MatView(mdn->A,viewer); CHKERRQ(ierr);
5128965ea79SLois Curfman McInnes     }
5138965ea79SLois Curfman McInnes     else {
5148965ea79SLois Curfman McInnes       /* assemble the entire matrix onto first processor. */
5158965ea79SLois Curfman McInnes       Mat          A;
51639ddd567SLois Curfman McInnes       int          M = mdn->M, N = mdn->N,m,row,i, nz, *cols;
51739ddd567SLois Curfman McInnes       Scalar       *vals;
51839ddd567SLois Curfman McInnes       Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data;
5198965ea79SLois Curfman McInnes 
5208965ea79SLois Curfman McInnes       if (!rank) {
521b4fd4287SBarry Smith         ierr = MatCreateMPIDense(mat->comm,M,M,N,N,PETSC_NULL,&A); CHKERRQ(ierr);
5228965ea79SLois Curfman McInnes       }
5238965ea79SLois Curfman McInnes       else {
524b4fd4287SBarry Smith         ierr = MatCreateMPIDense(mat->comm,0,M,N,N,PETSC_NULL,&A); CHKERRQ(ierr);
5258965ea79SLois Curfman McInnes       }
5268965ea79SLois Curfman McInnes       PLogObjectParent(mat,A);
5278965ea79SLois Curfman McInnes 
52839ddd567SLois Curfman McInnes       /* Copy the matrix ... This isn't the most efficient means,
52939ddd567SLois Curfman McInnes          but it's quick for now */
53039ddd567SLois Curfman McInnes       row = mdn->rstart; m = Amdn->m;
5318965ea79SLois Curfman McInnes       for ( i=0; i<m; i++ ) {
53239ddd567SLois Curfman McInnes         ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
53339ddd567SLois Curfman McInnes         ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr);
53439ddd567SLois Curfman McInnes         ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr);
53539ddd567SLois Curfman McInnes         row++;
5368965ea79SLois Curfman McInnes       }
5378965ea79SLois Curfman McInnes 
5386d4a8577SBarry Smith       ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5396d4a8577SBarry Smith       ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5408965ea79SLois Curfman McInnes       if (!rank) {
54139ddd567SLois Curfman McInnes         ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr);
5428965ea79SLois Curfman McInnes       }
5438965ea79SLois Curfman McInnes       ierr = MatDestroy(A); CHKERRQ(ierr);
5448965ea79SLois Curfman McInnes     }
5458965ea79SLois Curfman McInnes   }
5468965ea79SLois Curfman McInnes   return 0;
5478965ea79SLois Curfman McInnes }
5488965ea79SLois Curfman McInnes 
54939ddd567SLois Curfman McInnes static int MatView_MPIDense(PetscObject obj,Viewer viewer)
5508965ea79SLois Curfman McInnes {
5518965ea79SLois Curfman McInnes   Mat          mat = (Mat) obj;
55239ddd567SLois Curfman McInnes   int          ierr;
553bcd2baecSBarry Smith   ViewerType   vtype;
5548965ea79SLois Curfman McInnes 
555bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
556bcd2baecSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
55739ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
5588965ea79SLois Curfman McInnes   }
559bcd2baecSBarry Smith   else if (vtype == ASCII_FILES_VIEWER) {
56039ddd567SLois Curfman McInnes     ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr);
5618965ea79SLois Curfman McInnes   }
562bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
56339ddd567SLois Curfman McInnes     return MatView_MPIDense_Binary(mat,viewer);
5648965ea79SLois Curfman McInnes   }
5658965ea79SLois Curfman McInnes   return 0;
5668965ea79SLois Curfman McInnes }
5678965ea79SLois Curfman McInnes 
5684e220ebcSLois Curfman McInnes static int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
5698965ea79SLois Curfman McInnes {
5703501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
5713501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
5724e220ebcSLois Curfman McInnes   int          ierr;
5734e220ebcSLois Curfman McInnes   double       isend[5], irecv[5];
5748965ea79SLois Curfman McInnes 
5754e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
5764e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
5774e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
5784e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
5794e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
5804e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info); CHKERRQ(ierr);
5814e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
5824e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
5838965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
5844e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
5854e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
5864e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
5874e220ebcSLois Curfman McInnes     info->memory       = isend[3];
5884e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
5898965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
5903501a2bdSLois Curfman McInnes     MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_MAX,A->comm);
5914e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
5924e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
5934e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
5944e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
5954e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
5968965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
5973501a2bdSLois Curfman McInnes     MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_SUM,A->comm);
5984e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
5994e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6004e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6014e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6024e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6038965ea79SLois Curfman McInnes   }
6044e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6054e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6064e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6078965ea79SLois Curfman McInnes   return 0;
6088965ea79SLois Curfman McInnes }
6098965ea79SLois Curfman McInnes 
6108c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*);
6118aaee692SLois Curfman McInnes    extern int MatLUFactorNumeric_MPIDense(Mat,Mat*);
6128aaee692SLois Curfman McInnes    extern int MatLUFactor_MPIDense(Mat,IS,IS,double);
6138aaee692SLois Curfman McInnes    extern int MatSolve_MPIDense(Mat,Vec,Vec);
6148c469469SLois Curfman McInnes    extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec);
6158aaee692SLois Curfman McInnes    extern int MatSolveTrans_MPIDense(Mat,Vec,Vec);
6168aaee692SLois Curfman McInnes    extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */
6178aaee692SLois Curfman McInnes 
61839ddd567SLois Curfman McInnes static int MatSetOption_MPIDense(Mat A,MatOption op)
6198965ea79SLois Curfman McInnes {
62039ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
6218965ea79SLois Curfman McInnes 
6226d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
6236d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
6246d4a8577SBarry Smith       op == MAT_COLUMNS_SORTED ||
6256d4a8577SBarry Smith       op == MAT_ROW_ORIENTED) {
6268965ea79SLois Curfman McInnes         MatSetOption(a->A,op);
6278965ea79SLois Curfman McInnes   }
6286d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
6296d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
6306d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
6316d4a8577SBarry Smith            op == MAT_YES_NEW_DIAGONALS)
63294a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_MPIDense:Option ignored\n");
6336d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)
63439b7565bSBarry Smith     { a->roworiented = 0; MatSetOption(a->A,op);}
6356d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
6366d4a8577SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:MAT_NO_NEW_DIAGONALS");}
6378965ea79SLois Curfman McInnes   else
63839ddd567SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:unknown option");}
6398965ea79SLois Curfman McInnes   return 0;
6408965ea79SLois Curfman McInnes }
6418965ea79SLois Curfman McInnes 
6423501a2bdSLois Curfman McInnes static int MatGetSize_MPIDense(Mat A,int *m,int *n)
6438965ea79SLois Curfman McInnes {
6443501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6458965ea79SLois Curfman McInnes   *m = mat->M; *n = mat->N;
6468965ea79SLois Curfman McInnes   return 0;
6478965ea79SLois Curfman McInnes }
6488965ea79SLois Curfman McInnes 
6493501a2bdSLois Curfman McInnes static int MatGetLocalSize_MPIDense(Mat A,int *m,int *n)
6508965ea79SLois Curfman McInnes {
6513501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6528965ea79SLois Curfman McInnes   *m = mat->m; *n = mat->N;
6538965ea79SLois Curfman McInnes   return 0;
6548965ea79SLois Curfman McInnes }
6558965ea79SLois Curfman McInnes 
6563501a2bdSLois Curfman McInnes static int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n)
6578965ea79SLois Curfman McInnes {
6583501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
6598965ea79SLois Curfman McInnes   *m = mat->rstart; *n = mat->rend;
6608965ea79SLois Curfman McInnes   return 0;
6618965ea79SLois Curfman McInnes }
6628965ea79SLois Curfman McInnes 
6633501a2bdSLois Curfman McInnes static int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v)
6648965ea79SLois Curfman McInnes {
6653501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense *) A->data;
66639ddd567SLois Curfman McInnes   int          lrow, rstart = mat->rstart, rend = mat->rend;
6678965ea79SLois Curfman McInnes 
66839ddd567SLois Curfman McInnes   if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIDense:only local rows")
6698965ea79SLois Curfman McInnes   lrow = row - rstart;
67039ddd567SLois Curfman McInnes   return MatGetRow(mat->A,lrow,nz,idx,v);
6718965ea79SLois Curfman McInnes }
6728965ea79SLois Curfman McInnes 
67339ddd567SLois Curfman McInnes static int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v)
6748965ea79SLois Curfman McInnes {
6750452661fSBarry Smith   if (idx) PetscFree(*idx);
6760452661fSBarry Smith   if (v) PetscFree(*v);
6778965ea79SLois Curfman McInnes   return 0;
6788965ea79SLois Curfman McInnes }
6798965ea79SLois Curfman McInnes 
680cddf8d76SBarry Smith static int MatNorm_MPIDense(Mat A,NormType type,double *norm)
681096963f5SLois Curfman McInnes {
6823501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense *) A->data;
6833501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data;
6843501a2bdSLois Curfman McInnes   int          ierr, i, j;
6853501a2bdSLois Curfman McInnes   double       sum = 0.0;
6863501a2bdSLois Curfman McInnes   Scalar       *v = mat->v;
6873501a2bdSLois Curfman McInnes 
6883501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
6893501a2bdSLois Curfman McInnes     ierr =  MatNorm(mdn->A,type,norm); CHKERRQ(ierr);
6903501a2bdSLois Curfman McInnes   } else {
6913501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
6923501a2bdSLois Curfman McInnes       for (i=0; i<mat->n*mat->m; i++ ) {
6933501a2bdSLois Curfman McInnes #if defined(PETSC_COMPLEX)
6943501a2bdSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
6953501a2bdSLois Curfman McInnes #else
6963501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
6973501a2bdSLois Curfman McInnes #endif
6983501a2bdSLois Curfman McInnes       }
6996d84be18SBarry Smith       MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);
7003501a2bdSLois Curfman McInnes       *norm = sqrt(*norm);
7013501a2bdSLois Curfman McInnes       PLogFlops(2*mat->n*mat->m);
7023501a2bdSLois Curfman McInnes     }
7033501a2bdSLois Curfman McInnes     else if (type == NORM_1) {
7043501a2bdSLois Curfman McInnes       double *tmp, *tmp2;
7050452661fSBarry Smith       tmp  = (double *) PetscMalloc( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp);
7063501a2bdSLois Curfman McInnes       tmp2 = tmp + mdn->N;
707cddf8d76SBarry Smith       PetscMemzero(tmp,2*mdn->N*sizeof(double));
708096963f5SLois Curfman McInnes       *norm = 0.0;
7093501a2bdSLois Curfman McInnes       v = mat->v;
7103501a2bdSLois Curfman McInnes       for ( j=0; j<mat->n; j++ ) {
7113501a2bdSLois Curfman McInnes         for ( i=0; i<mat->m; i++ ) {
71267e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
7133501a2bdSLois Curfman McInnes         }
7143501a2bdSLois Curfman McInnes       }
7156d84be18SBarry Smith       MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);
7163501a2bdSLois Curfman McInnes       for ( j=0; j<mdn->N; j++ ) {
7173501a2bdSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
7183501a2bdSLois Curfman McInnes       }
7190452661fSBarry Smith       PetscFree(tmp);
7203501a2bdSLois Curfman McInnes       PLogFlops(mat->n*mat->m);
7213501a2bdSLois Curfman McInnes     }
7223501a2bdSLois Curfman McInnes     else if (type == NORM_INFINITY) { /* max row norm */
7233501a2bdSLois Curfman McInnes       double ntemp;
7243501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr);
7256d84be18SBarry Smith       MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);
7263501a2bdSLois Curfman McInnes     }
7273501a2bdSLois Curfman McInnes     else {
7283501a2bdSLois Curfman McInnes       SETERRQ(1,"MatNorm_MPIDense:No support for two norm");
7293501a2bdSLois Curfman McInnes     }
7303501a2bdSLois Curfman McInnes   }
7313501a2bdSLois Curfman McInnes   return 0;
7323501a2bdSLois Curfman McInnes }
7333501a2bdSLois Curfman McInnes 
7343501a2bdSLois Curfman McInnes static int MatTranspose_MPIDense(Mat A,Mat *matout)
7353501a2bdSLois Curfman McInnes {
7363501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense *) A->data;
7373501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data;
7383501a2bdSLois Curfman McInnes   Mat          B;
7393501a2bdSLois Curfman McInnes   int          M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart;
7403501a2bdSLois Curfman McInnes   int          j, i, ierr;
7413501a2bdSLois Curfman McInnes   Scalar       *v;
7423501a2bdSLois Curfman McInnes 
7437056b6fcSBarry Smith   if (matout == PETSC_NULL && M != N) {
7443501a2bdSLois Curfman McInnes     SETERRQ(1,"MatTranspose_MPIDense:Supports square matrix only in-place");
7457056b6fcSBarry Smith   }
7467056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
7473501a2bdSLois Curfman McInnes 
7483501a2bdSLois Curfman McInnes   m = Aloc->m; n = Aloc->n; v = Aloc->v;
7490452661fSBarry Smith   rwork = (int *) PetscMalloc(n*sizeof(int)); CHKPTRQ(rwork);
7503501a2bdSLois Curfman McInnes   for ( j=0; j<n; j++ ) {
7513501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
7523501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr);
7533501a2bdSLois Curfman McInnes     v   += m;
7543501a2bdSLois Curfman McInnes   }
7550452661fSBarry Smith   PetscFree(rwork);
7566d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
7576d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
7583638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
7593501a2bdSLois Curfman McInnes     *matout = B;
7603501a2bdSLois Curfman McInnes   } else {
7613501a2bdSLois Curfman McInnes     /* This isn't really an in-place transpose, but free data struct from a */
7620452661fSBarry Smith     PetscFree(a->rowners);
7633501a2bdSLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
7643501a2bdSLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
7653501a2bdSLois Curfman McInnes     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
7660452661fSBarry Smith     PetscFree(a);
7673501a2bdSLois Curfman McInnes     PetscMemcpy(A,B,sizeof(struct _Mat));
7680452661fSBarry Smith     PetscHeaderDestroy(B);
7693501a2bdSLois Curfman McInnes   }
770096963f5SLois Curfman McInnes   return 0;
771096963f5SLois Curfman McInnes }
772096963f5SLois Curfman McInnes 
77344cd7ae7SLois Curfman McInnes #include "pinclude/plapack.h"
77444cd7ae7SLois Curfman McInnes static int MatScale_MPIDense(Scalar *alpha,Mat inA)
77544cd7ae7SLois Curfman McInnes {
77644cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense *) inA->data;
77744cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->A->data;
77844cd7ae7SLois Curfman McInnes   int          one = 1, nz;
77944cd7ae7SLois Curfman McInnes 
78044cd7ae7SLois Curfman McInnes   nz = a->m*a->n;
78144cd7ae7SLois Curfman McInnes   BLscal_( &nz, alpha, a->v, &one );
78244cd7ae7SLois Curfman McInnes   PLogFlops(nz);
78344cd7ae7SLois Curfman McInnes   return 0;
78444cd7ae7SLois Curfman McInnes }
78544cd7ae7SLois Curfman McInnes 
7863d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat,Mat *,int);
7878965ea79SLois Curfman McInnes 
7888965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
78939ddd567SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIDense,
79039ddd567SLois Curfman McInnes        MatGetRow_MPIDense,MatRestoreRow_MPIDense,
79139ddd567SLois Curfman McInnes        MatMult_MPIDense,MatMultAdd_MPIDense,
792096963f5SLois Curfman McInnes        MatMultTrans_MPIDense,MatMultTransAdd_MPIDense,
793e7ca0642SLois Curfman McInnes /*       MatSolve_MPIDense,0, */
794e7ca0642SLois Curfman McInnes        0,0,
79539ddd567SLois Curfman McInnes        0,0,
7968c469469SLois Curfman McInnes        0,0,
7978c469469SLois Curfman McInnes /*       MatLUFactor_MPIDense,0, */
7988aaee692SLois Curfman McInnes        0,MatTranspose_MPIDense,
79939ddd567SLois Curfman McInnes        MatGetInfo_MPIDense,0,
800096963f5SLois Curfman McInnes        MatGetDiagonal_MPIDense,0,MatNorm_MPIDense,
80139ddd567SLois Curfman McInnes        MatAssemblyBegin_MPIDense,MatAssemblyEnd_MPIDense,
8028965ea79SLois Curfman McInnes        0,
80339ddd567SLois Curfman McInnes        MatSetOption_MPIDense,MatZeroEntries_MPIDense,MatZeroRows_MPIDense,
8043b2fbd54SBarry Smith        0,0,
8058c469469SLois Curfman McInnes /*       0,MatLUFactorSymbolic_MPIDense,MatLUFactorNumeric_MPIDense, */
8068aaee692SLois Curfman McInnes        0,0,
80739ddd567SLois Curfman McInnes        MatGetSize_MPIDense,MatGetLocalSize_MPIDense,
80839ddd567SLois Curfman McInnes        MatGetOwnershipRange_MPIDense,
809ff14e315SSatish Balay        0,0, MatGetArray_MPIDense, MatRestoreArray_MPIDense,
810905e6a2fSBarry Smith        0,MatConvertSameType_MPIDense,
811b49de8d1SLois Curfman McInnes        0,0,0,0,0,
8124e220ebcSLois Curfman McInnes        0,0,MatGetValues_MPIDense,0,0,MatScale_MPIDense,
8134e220ebcSLois Curfman McInnes        0,0,0,MatGetBlockSize_MPIDense};
8148965ea79SLois Curfman McInnes 
8158965ea79SLois Curfman McInnes /*@C
81639ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
8178965ea79SLois Curfman McInnes 
8188965ea79SLois Curfman McInnes    Input Parameters:
8198965ea79SLois Curfman McInnes .  comm - MPI communicator
8208965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
8218965ea79SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated
8228965ea79SLois Curfman McInnes            if N is given)
8238965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
8248965ea79SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated
8258965ea79SLois Curfman McInnes            if n is given)
826b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
827dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
8288965ea79SLois Curfman McInnes 
8298965ea79SLois Curfman McInnes    Output Parameter:
830477f1c0bSLois Curfman McInnes .  A - the matrix
8318965ea79SLois Curfman McInnes 
8328965ea79SLois Curfman McInnes    Notes:
83339ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
83439ddd567SLois Curfman McInnes    storage by columns.
8358965ea79SLois Curfman McInnes 
83618f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
83718f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
838b4fd4287SBarry Smith    set data=PETSC_NULL.
83918f449edSLois Curfman McInnes 
8408965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
8418965ea79SLois Curfman McInnes    (possibly both).
8428965ea79SLois Curfman McInnes 
8433501a2bdSLois Curfman McInnes    Currently, the only parallel dense matrix decomposition is by rows,
8443501a2bdSLois Curfman McInnes    so that n=N and each submatrix owns all of the global columns.
8453501a2bdSLois Curfman McInnes 
84639ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel
8478965ea79SLois Curfman McInnes 
84839ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
8498965ea79SLois Curfman McInnes @*/
850477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A)
8518965ea79SLois Curfman McInnes {
8528965ea79SLois Curfman McInnes   Mat          mat;
85339ddd567SLois Curfman McInnes   Mat_MPIDense *a;
85425cdf11fSBarry Smith   int          ierr, i,flg;
8558965ea79SLois Curfman McInnes 
856ed2daf61SLois Curfman McInnes   /* Note:  For now, when data is specified above, this assumes the user correctly
857ed2daf61SLois Curfman McInnes    allocates the local dense storage space.  We should add error checking. */
85818f449edSLois Curfman McInnes 
859477f1c0bSLois Curfman McInnes   *A = 0;
8600452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIDENSE,comm);
8618965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
8620452661fSBarry Smith   mat->data       = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
8638965ea79SLois Curfman McInnes   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
86439ddd567SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIDense;
86539ddd567SLois Curfman McInnes   mat->view       = MatView_MPIDense;
8668965ea79SLois Curfman McInnes   mat->factor     = 0;
8678965ea79SLois Curfman McInnes 
868622d7880SLois Curfman McInnes   a->factor     = 0;
8698965ea79SLois Curfman McInnes   a->insertmode = NOT_SET_VALUES;
8708965ea79SLois Curfman McInnes   MPI_Comm_rank(comm,&a->rank);
8718965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&a->size);
8728965ea79SLois Curfman McInnes 
87339ddd567SLois Curfman McInnes   if (M == PETSC_DECIDE) MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm);
8748965ea79SLois Curfman McInnes   if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);}
87539ddd567SLois Curfman McInnes 
87639ddd567SLois Curfman McInnes   /* each row stores all columns */
87739ddd567SLois Curfman McInnes   if (N == PETSC_DECIDE) N = n;
878d7e8b826SBarry Smith   if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);}
879d7e8b826SBarry Smith   /*  if (n != N) SETERRQ(1,"MatCreateMPIDense:For now, only n=N is supported"); */
880aca0ad90SLois Curfman McInnes   a->N = mat->N = N;
881aca0ad90SLois Curfman McInnes   a->M = mat->M = M;
882aca0ad90SLois Curfman McInnes   a->m = mat->m = m;
883aca0ad90SLois Curfman McInnes   a->n = mat->n = n;
8848965ea79SLois Curfman McInnes 
8858965ea79SLois Curfman McInnes   /* build local table of row and column ownerships */
886d7e8b826SBarry Smith   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners);
887d7e8b826SBarry Smith   a->cowners = a->rowners + a->size + 1;
8887056b6fcSBarry Smith   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIDense));
8898965ea79SLois Curfman McInnes   MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);
8908965ea79SLois Curfman McInnes   a->rowners[0] = 0;
8918965ea79SLois Curfman McInnes   for ( i=2; i<=a->size; i++ ) {
8928965ea79SLois Curfman McInnes     a->rowners[i] += a->rowners[i-1];
8938965ea79SLois Curfman McInnes   }
8948965ea79SLois Curfman McInnes   a->rstart = a->rowners[a->rank];
8958965ea79SLois Curfman McInnes   a->rend   = a->rowners[a->rank+1];
896d7e8b826SBarry Smith   MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);
897d7e8b826SBarry Smith   a->cowners[0] = 0;
898d7e8b826SBarry Smith   for ( i=2; i<=a->size; i++ ) {
899d7e8b826SBarry Smith     a->cowners[i] += a->cowners[i-1];
900d7e8b826SBarry Smith   }
9018965ea79SLois Curfman McInnes 
90218f449edSLois Curfman McInnes   ierr = MatCreateSeqDense(MPI_COMM_SELF,m,N,data,&a->A); CHKERRQ(ierr);
9038965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
9048965ea79SLois Curfman McInnes 
9058965ea79SLois Curfman McInnes   /* build cache for off array entries formed */
9068965ea79SLois Curfman McInnes   ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr);
9078965ea79SLois Curfman McInnes 
9088965ea79SLois Curfman McInnes   /* stuff used for matrix vector multiply */
9098965ea79SLois Curfman McInnes   a->lvec        = 0;
9108965ea79SLois Curfman McInnes   a->Mvctx       = 0;
91139b7565bSBarry Smith   a->roworiented = 1;
9128965ea79SLois Curfman McInnes 
913477f1c0bSLois Curfman McInnes   *A = mat;
91425cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
91525cdf11fSBarry Smith   if (flg) {
9168c469469SLois Curfman McInnes     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
9178c469469SLois Curfman McInnes   }
9188965ea79SLois Curfman McInnes   return 0;
9198965ea79SLois Curfman McInnes }
9208965ea79SLois Curfman McInnes 
9213d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat A,Mat *newmat,int cpvalues)
9228965ea79SLois Curfman McInnes {
9238965ea79SLois Curfman McInnes   Mat          mat;
9243501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data;
92539ddd567SLois Curfman McInnes   int          ierr;
9262ba99913SLois Curfman McInnes   FactorCtx    *factor;
9278965ea79SLois Curfman McInnes 
9288965ea79SLois Curfman McInnes   *newmat       = 0;
9290452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIDENSE,A->comm);
9308965ea79SLois Curfman McInnes   PLogObjectCreate(mat);
9310452661fSBarry Smith   mat->data      = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a);
9328965ea79SLois Curfman McInnes   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
93339ddd567SLois Curfman McInnes   mat->destroy   = MatDestroy_MPIDense;
93439ddd567SLois Curfman McInnes   mat->view      = MatView_MPIDense;
9353501a2bdSLois Curfman McInnes   mat->factor    = A->factor;
936c456f294SBarry Smith   mat->assembled = PETSC_TRUE;
9378965ea79SLois Curfman McInnes 
93844cd7ae7SLois Curfman McInnes   a->m = mat->m = oldmat->m;
93944cd7ae7SLois Curfman McInnes   a->n = mat->n = oldmat->n;
94044cd7ae7SLois Curfman McInnes   a->M = mat->M = oldmat->M;
94144cd7ae7SLois Curfman McInnes   a->N = mat->N = oldmat->N;
9422ba99913SLois Curfman McInnes   if (oldmat->factor) {
9432ba99913SLois Curfman McInnes     a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx)); CHKPTRQ(factor);
9442ba99913SLois Curfman McInnes     /* copy factor contents ... add this code! */
9452ba99913SLois Curfman McInnes   } else a->factor = 0;
9468965ea79SLois Curfman McInnes 
9478965ea79SLois Curfman McInnes   a->rstart     = oldmat->rstart;
9488965ea79SLois Curfman McInnes   a->rend       = oldmat->rend;
9498965ea79SLois Curfman McInnes   a->size       = oldmat->size;
9508965ea79SLois Curfman McInnes   a->rank       = oldmat->rank;
9518965ea79SLois Curfman McInnes   a->insertmode = NOT_SET_VALUES;
9528965ea79SLois Curfman McInnes 
9530452661fSBarry Smith   a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners);
95439ddd567SLois Curfman McInnes   PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIDense));
9558965ea79SLois Curfman McInnes   PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));
9568965ea79SLois Curfman McInnes   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
9578965ea79SLois Curfman McInnes 
9588965ea79SLois Curfman McInnes   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
9598965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->lvec);
96055659b69SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
9618965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->Mvctx);
9628965ea79SLois Curfman McInnes   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
9638965ea79SLois Curfman McInnes   PLogObjectParent(mat,a->A);
9648965ea79SLois Curfman McInnes   *newmat = mat;
9658965ea79SLois Curfman McInnes   return 0;
9668965ea79SLois Curfman McInnes }
9678965ea79SLois Curfman McInnes 
96877c4ece6SBarry Smith #include "sys.h"
9698965ea79SLois Curfman McInnes 
97090ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat)
97190ace30eSBarry Smith {
97290ace30eSBarry Smith   int        *rowners, i,size,rank,m,rstart,rend,ierr,nz,j;
97390ace30eSBarry Smith   Scalar     *array,*vals,*vals_ptr;
97490ace30eSBarry Smith   MPI_Status status;
97590ace30eSBarry Smith 
97690ace30eSBarry Smith   MPI_Comm_rank(comm,&rank);
97790ace30eSBarry Smith   MPI_Comm_size(comm,&size);
97890ace30eSBarry Smith 
97990ace30eSBarry Smith   /* determine ownership of all rows */
98090ace30eSBarry Smith   m = M/size + ((M % size) > rank);
98190ace30eSBarry Smith   rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
98290ace30eSBarry Smith   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
98390ace30eSBarry Smith   rowners[0] = 0;
98490ace30eSBarry Smith   for ( i=2; i<=size; i++ ) {
98590ace30eSBarry Smith     rowners[i] += rowners[i-1];
98690ace30eSBarry Smith   }
98790ace30eSBarry Smith   rstart = rowners[rank];
98890ace30eSBarry Smith   rend   = rowners[rank+1];
98990ace30eSBarry Smith 
99090ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
99190ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array); CHKERRQ(ierr);
99290ace30eSBarry Smith 
99390ace30eSBarry Smith   if (!rank) {
99490ace30eSBarry Smith     vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
99590ace30eSBarry Smith 
99690ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
99777c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,BINARY_SCALAR); CHKERRQ(ierr);
99890ace30eSBarry Smith 
99990ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
100090ace30eSBarry Smith     vals_ptr = vals;
100190ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
100290ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
100390ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
100490ace30eSBarry Smith       }
100590ace30eSBarry Smith     }
100690ace30eSBarry Smith 
100790ace30eSBarry Smith     /* read in other processors and ship out */
100890ace30eSBarry Smith     for ( i=1; i<size; i++ ) {
100990ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
101077c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
101190ace30eSBarry Smith       MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);
101290ace30eSBarry Smith     }
101390ace30eSBarry Smith   }
101490ace30eSBarry Smith   else {
101590ace30eSBarry Smith     /* receive numeric values */
101690ace30eSBarry Smith     vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals);
101790ace30eSBarry Smith 
101890ace30eSBarry Smith     /* receive message of values*/
101990ace30eSBarry Smith     MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);
102090ace30eSBarry Smith 
102190ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
102290ace30eSBarry Smith     vals_ptr = vals;
102390ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
102490ace30eSBarry Smith       for ( j=0; j<N; j++ ) {
102590ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
102690ace30eSBarry Smith       }
102790ace30eSBarry Smith     }
102890ace30eSBarry Smith   }
102990ace30eSBarry Smith   PetscFree(rowners);
103090ace30eSBarry Smith   PetscFree(vals);
10316d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10326d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
103390ace30eSBarry Smith   return 0;
103490ace30eSBarry Smith }
103590ace30eSBarry Smith 
103690ace30eSBarry Smith 
103719bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat)
10388965ea79SLois Curfman McInnes {
10398965ea79SLois Curfman McInnes   Mat          A;
10408965ea79SLois Curfman McInnes   int          i, nz, ierr, j,rstart, rend, fd;
10418965ea79SLois Curfman McInnes   Scalar       *vals,*svals;
104219bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
10438965ea79SLois Curfman McInnes   MPI_Status   status;
10448965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
10458965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
104619bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
10478965ea79SLois Curfman McInnes 
10488965ea79SLois Curfman McInnes   MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank);
10498965ea79SLois Curfman McInnes   if (!rank) {
105090ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
105177c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr);
105239ddd567SLois Curfman McInnes     if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:not matrix object");
10538965ea79SLois Curfman McInnes   }
10548965ea79SLois Curfman McInnes 
10558965ea79SLois Curfman McInnes   MPI_Bcast(header+1,3,MPI_INT,0,comm);
105690ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
105790ace30eSBarry Smith 
105890ace30eSBarry Smith   /*
105990ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
106090ace30eSBarry Smith   */
106190ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
106290ace30eSBarry Smith     return MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);
106390ace30eSBarry Smith   }
106490ace30eSBarry Smith 
10658965ea79SLois Curfman McInnes   /* determine ownership of all rows */
10668965ea79SLois Curfman McInnes   m = M/size + ((M % size) > rank);
10670452661fSBarry Smith   rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
10688965ea79SLois Curfman McInnes   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
10698965ea79SLois Curfman McInnes   rowners[0] = 0;
10708965ea79SLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
10718965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
10728965ea79SLois Curfman McInnes   }
10738965ea79SLois Curfman McInnes   rstart = rowners[rank];
10748965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
10758965ea79SLois Curfman McInnes 
10768965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
10770452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
10788965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
10798965ea79SLois Curfman McInnes   if (!rank) {
10800452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
108177c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
10820452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts);
10838965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
10848965ea79SLois Curfman McInnes     MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
10850452661fSBarry Smith     PetscFree(sndcounts);
10868965ea79SLois Curfman McInnes   }
10878965ea79SLois Curfman McInnes   else {
10888965ea79SLois Curfman McInnes     MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);
10898965ea79SLois Curfman McInnes   }
10908965ea79SLois Curfman McInnes 
10918965ea79SLois Curfman McInnes   if (!rank) {
10928965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
10930452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz);
1094cddf8d76SBarry Smith     PetscMemzero(procsnz,size*sizeof(int));
10958965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
10968965ea79SLois Curfman McInnes       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
10978965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
10988965ea79SLois Curfman McInnes       }
10998965ea79SLois Curfman McInnes     }
11000452661fSBarry Smith     PetscFree(rowlengths);
11018965ea79SLois Curfman McInnes 
11028965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
11038965ea79SLois Curfman McInnes     maxnz = 0;
11048965ea79SLois Curfman McInnes     for ( i=0; i<size; i++ ) {
11050452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
11068965ea79SLois Curfman McInnes     }
11070452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols);
11088965ea79SLois Curfman McInnes 
11098965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
11108965ea79SLois Curfman McInnes     nz = procsnz[0];
11110452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
111277c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr);
11138965ea79SLois Curfman McInnes 
11148965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
11158965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
11168965ea79SLois Curfman McInnes       nz = procsnz[i];
111777c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
11188965ea79SLois Curfman McInnes       MPI_Send(cols,nz,MPI_INT,i,tag,comm);
11198965ea79SLois Curfman McInnes     }
11200452661fSBarry Smith     PetscFree(cols);
11218965ea79SLois Curfman McInnes   }
11228965ea79SLois Curfman McInnes   else {
11238965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
11248965ea79SLois Curfman McInnes     nz = 0;
11258965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
11268965ea79SLois Curfman McInnes       nz += ourlens[i];
11278965ea79SLois Curfman McInnes     }
11280452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
11298965ea79SLois Curfman McInnes 
11308965ea79SLois Curfman McInnes     /* receive message of column indices*/
11318965ea79SLois Curfman McInnes     MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);
11328965ea79SLois Curfman McInnes     MPI_Get_count(&status,MPI_INT,&maxnz);
113339ddd567SLois Curfman McInnes     if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file");
11348965ea79SLois Curfman McInnes   }
11358965ea79SLois Curfman McInnes 
11368965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1137cddf8d76SBarry Smith   PetscMemzero(offlens,m*sizeof(int));
11388965ea79SLois Curfman McInnes   jj = 0;
11398965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
11408965ea79SLois Curfman McInnes     for ( j=0; j<ourlens[i]; j++ ) {
11418965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
11428965ea79SLois Curfman McInnes       jj++;
11438965ea79SLois Curfman McInnes     }
11448965ea79SLois Curfman McInnes   }
11458965ea79SLois Curfman McInnes 
11468965ea79SLois Curfman McInnes   /* create our matrix */
11478965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
11488965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
11498965ea79SLois Curfman McInnes   }
1150b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
11518965ea79SLois Curfman McInnes   A = *newmat;
11528965ea79SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
11538965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
11548965ea79SLois Curfman McInnes   }
11558965ea79SLois Curfman McInnes 
11568965ea79SLois Curfman McInnes   if (!rank) {
11570452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
11588965ea79SLois Curfman McInnes 
11598965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
11608965ea79SLois Curfman McInnes     nz = procsnz[0];
116177c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
11628965ea79SLois Curfman McInnes 
11638965ea79SLois Curfman McInnes     /* insert into matrix */
11648965ea79SLois Curfman McInnes     jj      = rstart;
11658965ea79SLois Curfman McInnes     smycols = mycols;
11668965ea79SLois Curfman McInnes     svals   = vals;
11678965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
11688965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
11698965ea79SLois Curfman McInnes       smycols += ourlens[i];
11708965ea79SLois Curfman McInnes       svals   += ourlens[i];
11718965ea79SLois Curfman McInnes       jj++;
11728965ea79SLois Curfman McInnes     }
11738965ea79SLois Curfman McInnes 
11748965ea79SLois Curfman McInnes     /* read in other processors and ship out */
11758965ea79SLois Curfman McInnes     for ( i=1; i<size; i++ ) {
11768965ea79SLois Curfman McInnes       nz = procsnz[i];
117777c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
11788965ea79SLois Curfman McInnes       MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
11798965ea79SLois Curfman McInnes     }
11800452661fSBarry Smith     PetscFree(procsnz);
11818965ea79SLois Curfman McInnes   }
11828965ea79SLois Curfman McInnes   else {
11838965ea79SLois Curfman McInnes     /* receive numeric values */
11840452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
11858965ea79SLois Curfman McInnes 
11868965ea79SLois Curfman McInnes     /* receive message of values*/
11878965ea79SLois Curfman McInnes     MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
11888965ea79SLois Curfman McInnes     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
118939ddd567SLois Curfman McInnes     if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file");
11908965ea79SLois Curfman McInnes 
11918965ea79SLois Curfman McInnes     /* insert into matrix */
11928965ea79SLois Curfman McInnes     jj      = rstart;
11938965ea79SLois Curfman McInnes     smycols = mycols;
11948965ea79SLois Curfman McInnes     svals   = vals;
11958965ea79SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
11968965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
11978965ea79SLois Curfman McInnes       smycols += ourlens[i];
11988965ea79SLois Curfman McInnes       svals   += ourlens[i];
11998965ea79SLois Curfman McInnes       jj++;
12008965ea79SLois Curfman McInnes     }
12018965ea79SLois Curfman McInnes   }
12020452661fSBarry Smith   PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners);
12038965ea79SLois Curfman McInnes 
12046d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12056d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12068965ea79SLois Curfman McInnes   return 0;
12078965ea79SLois Curfman McInnes }
120890ace30eSBarry Smith 
120990ace30eSBarry Smith 
121090ace30eSBarry Smith 
121190ace30eSBarry Smith 
121290ace30eSBarry Smith 
1213