xref: /petsc/src/mat/impls/aij/mpi/mpiaij.c (revision 3369ce9a255b285539dc954b9ad9cefcea2f84b6)
1cb512458SBarry Smith #ifndef lint
2*3369ce9aSBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.194 1997/03/13 16:33:58 curfman Exp bsmith $";
3cb512458SBarry Smith #endif
48a729477SBarry Smith 
5*3369ce9aSBarry Smith #include "pinclude/pviewer.h"
670f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h"
7f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
8d9942c19SSatish Balay #include "src/inline/spops.h"
98a729477SBarry Smith 
109e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column
119e25ed09SBarry Smith number to the local number in the off-diagonal part of the local
129e25ed09SBarry Smith storage of the matrix.  This is done in a non scable way since the
139e25ed09SBarry Smith length of colmap equals the global matrix length.
149e25ed09SBarry Smith */
155615d1e5SSatish Balay #undef __FUNC__
165eea60f9SBarry Smith #define __FUNC__ "CreateColmap_MPIAIJ_Private" /* ADIC Ignore */
170a198c4cSBarry Smith int CreateColmap_MPIAIJ_Private(Mat mat)
189e25ed09SBarry Smith {
1944a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
20ec8511deSBarry Smith   Mat_SeqAIJ *B = (Mat_SeqAIJ*) aij->B->data;
21905e6a2fSBarry Smith   int        n = B->n,i;
22dbb450caSBarry Smith 
230452661fSBarry Smith   aij->colmap = (int *) PetscMalloc(aij->N*sizeof(int));CHKPTRQ(aij->colmap);
24464493b3SBarry Smith   PLogObjectMemory(mat,aij->N*sizeof(int));
25cddf8d76SBarry Smith   PetscMemzero(aij->colmap,aij->N*sizeof(int));
26905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1;
279e25ed09SBarry Smith   return 0;
289e25ed09SBarry Smith }
299e25ed09SBarry Smith 
302493cbb0SBarry Smith extern int DisAssemble_MPIAIJ(Mat);
312493cbb0SBarry Smith 
325615d1e5SSatish Balay #undef __FUNC__
335615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_MPIAIJ"
348f6be9afSLois Curfman McInnes int MatGetRowIJ_MPIAIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int **ja,
353b2fbd54SBarry Smith                            PetscTruth *done)
36299609e3SLois Curfman McInnes {
37299609e3SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
38299609e3SLois Curfman McInnes   int        ierr;
3917699dbbSLois Curfman McInnes   if (aij->size == 1) {
403b2fbd54SBarry Smith     ierr = MatGetRowIJ(aij->A,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
41e3372554SBarry Smith   } else SETERRQ(1,0,"not supported in parallel");
423b2fbd54SBarry Smith   return 0;
433b2fbd54SBarry Smith }
443b2fbd54SBarry Smith 
455615d1e5SSatish Balay #undef __FUNC__
465eea60f9SBarry Smith #define __FUNC__ "MatRestoreRowIJ_MPIAIJ" /* ADIC Ignore */
478f6be9afSLois Curfman McInnes int MatRestoreRowIJ_MPIAIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int **ja,
483b2fbd54SBarry Smith                                PetscTruth *done)
493b2fbd54SBarry Smith {
503b2fbd54SBarry Smith   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
513b2fbd54SBarry Smith   int        ierr;
523b2fbd54SBarry Smith   if (aij->size == 1) {
533b2fbd54SBarry Smith     ierr = MatRestoreRowIJ(aij->A,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
54e3372554SBarry Smith   } else SETERRQ(1,0,"not supported in parallel");
55299609e3SLois Curfman McInnes   return 0;
56299609e3SLois Curfman McInnes }
57299609e3SLois Curfman McInnes 
580520107fSSatish Balay #define CHUNKSIZE   15
59f5e9677aSSatish Balay #define MatSetValues_SeqAIJ_A_Private(A,row,col,value,addv) \
600520107fSSatish Balay { \
610520107fSSatish Balay  \
620520107fSSatish Balay     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift; \
630520107fSSatish Balay     rmax = imax[row]; nrow = ailen[row];  \
64f5e9677aSSatish Balay     col1 = col - shift; \
65f5e9677aSSatish Balay      \
660520107fSSatish Balay       for ( _i=0; _i<nrow; _i++ ) { \
67f5e9677aSSatish Balay         if (rp[_i] > col1) break; \
68f5e9677aSSatish Balay         if (rp[_i] == col1) { \
690520107fSSatish Balay           if (addv == ADD_VALUES) ap[_i] += value;   \
700520107fSSatish Balay           else                  ap[_i] = value; \
710520107fSSatish Balay           goto _noinsert; \
720520107fSSatish Balay         } \
730520107fSSatish Balay       }  \
74f5e9677aSSatish Balay       if (nonew)  goto _noinsert; \
750520107fSSatish Balay       if (nrow >= rmax) { \
760520107fSSatish Balay         /* there is no extra room in row, therefore enlarge */ \
770520107fSSatish Balay         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; \
780520107fSSatish Balay         Scalar *new_a; \
790520107fSSatish Balay  \
800520107fSSatish Balay         /* malloc new storage space */ \
810520107fSSatish Balay         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); \
820520107fSSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \
830520107fSSatish Balay         new_j   = (int *) (new_a + new_nz); \
840520107fSSatish Balay         new_i   = new_j + new_nz; \
850520107fSSatish Balay  \
860520107fSSatish Balay         /* copy over old data into new slots */ \
870520107fSSatish Balay         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} \
880520107fSSatish Balay         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} \
890520107fSSatish Balay         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); \
900520107fSSatish Balay         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); \
910520107fSSatish Balay         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, \
920520107fSSatish Balay                                                            len*sizeof(int)); \
930520107fSSatish Balay         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); \
940520107fSSatish Balay         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, \
950520107fSSatish Balay                                                            len*sizeof(Scalar));  \
960520107fSSatish Balay         /* free up old matrix storage */ \
97f5e9677aSSatish Balay  \
980520107fSSatish Balay         PetscFree(a->a);  \
990520107fSSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} \
1000520107fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;  \
1010520107fSSatish Balay         a->singlemalloc = 1; \
1020520107fSSatish Balay  \
1030520107fSSatish Balay         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift; \
1040520107fSSatish Balay         rmax = imax[row] = imax[row] + CHUNKSIZE; \
1050520107fSSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \
1060520107fSSatish Balay         a->maxnz += CHUNKSIZE; \
1070520107fSSatish Balay         a->reallocs++; \
1080520107fSSatish Balay       } \
1090520107fSSatish Balay       N = nrow++ - 1; a->nz++; \
1100520107fSSatish Balay       /* shift up all the later entries in this row */ \
1110520107fSSatish Balay       for ( ii=N; ii>=_i; ii-- ) { \
1120520107fSSatish Balay         rp[ii+1] = rp[ii]; \
1130520107fSSatish Balay         ap[ii+1] = ap[ii]; \
1140520107fSSatish Balay       } \
115f5e9677aSSatish Balay       rp[_i] = col1;  \
1160520107fSSatish Balay       ap[_i] = value;  \
1170520107fSSatish Balay       _noinsert: ; \
1180520107fSSatish Balay       ailen[row] = nrow; \
1190520107fSSatish Balay }
1200a198c4cSBarry Smith 
1210520107fSSatish Balay extern int MatSetValues_SeqAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode);
1225615d1e5SSatish Balay #undef __FUNC__
1235615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIAIJ"
1248f6be9afSLois Curfman McInnes int MatSetValues_MPIAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv)
1258a729477SBarry Smith {
12644a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1274b0e389bSBarry Smith   Scalar     value;
1281eb62cbbSBarry Smith   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
1291eb62cbbSBarry Smith   int        cstart = aij->cstart, cend = aij->cend,row,col;
130905e6a2fSBarry Smith   int        roworiented = aij->roworiented;
1318a729477SBarry Smith 
1320520107fSSatish Balay   /* Some Variables required in the macro */
1334ee7247eSSatish Balay   Mat        A = aij->A;
1344ee7247eSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1354ee7247eSSatish Balay   int        *rp,ii,nrow,_i,rmax, N, col1;
1360520107fSSatish Balay   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen;
1374ee7247eSSatish Balay   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
1380520107fSSatish Balay   Scalar     *ap, *aa = a->a;
1394ee7247eSSatish Balay 
1408a729477SBarry Smith   for ( i=0; i<m; i++ ) {
1410a198c4cSBarry Smith #if defined(PETSC_BOPT_g)
142e3372554SBarry Smith     if (im[i] < 0) SETERRQ(1,0,"Negative row");
143e3372554SBarry Smith     if (im[i] >= aij->M) SETERRQ(1,0,"Row too large");
1440a198c4cSBarry Smith #endif
1454b0e389bSBarry Smith     if (im[i] >= rstart && im[i] < rend) {
1464b0e389bSBarry Smith       row = im[i] - rstart;
1471eb62cbbSBarry Smith       for ( j=0; j<n; j++ ) {
1484b0e389bSBarry Smith         if (in[j] >= cstart && in[j] < cend){
1494b0e389bSBarry Smith           col = in[j] - cstart;
1504b0e389bSBarry Smith           if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
151f5e9677aSSatish Balay           MatSetValues_SeqAIJ_A_Private(aij->A,row,col,value,addv);
1520520107fSSatish Balay           /* ierr = MatSetValues_SeqAIJ(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */
1531eb62cbbSBarry Smith         }
1540a198c4cSBarry Smith #if defined(PETSC_BOPT_g)
155e3372554SBarry Smith         else if (in[j] < 0) {SETERRQ(1,0,"Negative column");}
156e3372554SBarry Smith         else if (in[j] >= aij->N) {SETERRQ(1,0,"Col too large");}
1570a198c4cSBarry Smith #endif
1581eb62cbbSBarry Smith         else {
159227d817aSBarry Smith           if (mat->was_assembled) {
160905e6a2fSBarry Smith             if (!aij->colmap) {
161905e6a2fSBarry Smith               ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
162905e6a2fSBarry Smith             }
163905e6a2fSBarry Smith             col = aij->colmap[in[j]] - 1;
164ec8511deSBarry Smith             if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) {
1652493cbb0SBarry Smith               ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
1664b0e389bSBarry Smith               col =  in[j];
167d6dfbf8fSBarry Smith             }
1689e25ed09SBarry Smith           }
1694b0e389bSBarry Smith           else col = in[j];
1704b0e389bSBarry Smith           if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
1710a198c4cSBarry Smith           ierr = MatSetValues_SeqAIJ(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr);
1721eb62cbbSBarry Smith         }
1731eb62cbbSBarry Smith       }
1741eb62cbbSBarry Smith     }
1751eb62cbbSBarry Smith     else {
17690f02eecSBarry Smith       if (roworiented && !aij->donotstash) {
1774b0e389bSBarry Smith         ierr = StashValues_Private(&aij->stash,im[i],n,in,v+i*n,addv);CHKERRQ(ierr);
1784b0e389bSBarry Smith       }
1794b0e389bSBarry Smith       else {
18090f02eecSBarry Smith         if (!aij->donotstash) {
1814b0e389bSBarry Smith           row = im[i];
1824b0e389bSBarry Smith           for ( j=0; j<n; j++ ) {
1834b0e389bSBarry Smith             ierr = StashValues_Private(&aij->stash,row,1,in+j,v+i+j*m,addv);CHKERRQ(ierr);
1844b0e389bSBarry Smith           }
1854b0e389bSBarry Smith         }
1861eb62cbbSBarry Smith       }
1878a729477SBarry Smith     }
18890f02eecSBarry Smith   }
1898a729477SBarry Smith   return 0;
1908a729477SBarry Smith }
1918a729477SBarry Smith 
1925615d1e5SSatish Balay #undef __FUNC__
1935615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIAIJ"
1948f6be9afSLois Curfman McInnes int MatGetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
195b49de8d1SLois Curfman McInnes {
196b49de8d1SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
197b49de8d1SLois Curfman McInnes   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
198b49de8d1SLois Curfman McInnes   int        cstart = aij->cstart, cend = aij->cend,row,col;
199b49de8d1SLois Curfman McInnes 
200b49de8d1SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
201e3372554SBarry Smith     if (idxm[i] < 0) SETERRQ(1,0,"Negative row");
202e3372554SBarry Smith     if (idxm[i] >= aij->M) SETERRQ(1,0,"Row too large");
203b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
204b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
205b49de8d1SLois Curfman McInnes       for ( j=0; j<n; j++ ) {
206e3372554SBarry Smith         if (idxn[j] < 0) SETERRQ(1,0,"Negative column");
207e3372554SBarry Smith         if (idxn[j] >= aij->N) SETERRQ(1,0,"Col too large");
208b49de8d1SLois Curfman McInnes         if (idxn[j] >= cstart && idxn[j] < cend){
209b49de8d1SLois Curfman McInnes           col = idxn[j] - cstart;
210b49de8d1SLois Curfman McInnes           ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr);
211b49de8d1SLois Curfman McInnes         }
212b49de8d1SLois Curfman McInnes         else {
213905e6a2fSBarry Smith           if (!aij->colmap) {
214905e6a2fSBarry Smith             ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
215905e6a2fSBarry Smith           }
216905e6a2fSBarry Smith           col = aij->colmap[idxn[j]] - 1;
217e60e1c95SSatish Balay           if ((col < 0) || (aij->garray[col] != idxn[j])) *(v+i*n+j) = 0.0;
218d9d09a02SSatish Balay           else {
219b49de8d1SLois Curfman McInnes             ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr);
220b49de8d1SLois Curfman McInnes           }
221b49de8d1SLois Curfman McInnes         }
222b49de8d1SLois Curfman McInnes       }
223d9d09a02SSatish Balay     }
224b49de8d1SLois Curfman McInnes     else {
225e3372554SBarry Smith       SETERRQ(1,0,"Only local values currently supported");
226b49de8d1SLois Curfman McInnes     }
227b49de8d1SLois Curfman McInnes   }
228b49de8d1SLois Curfman McInnes   return 0;
229b49de8d1SLois Curfman McInnes }
230b49de8d1SLois Curfman McInnes 
2315615d1e5SSatish Balay #undef __FUNC__
2325615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ"
2338f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode)
2348a729477SBarry Smith {
23544a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
236d6dfbf8fSBarry Smith   MPI_Comm    comm = mat->comm;
23717699dbbSLois Curfman McInnes   int         size = aij->size, *owners = aij->rowners;
23817699dbbSLois Curfman McInnes   int         rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr;
2391eb62cbbSBarry Smith   MPI_Request *send_waits,*recv_waits;
2406abc6512SBarry Smith   int         *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
2411eb62cbbSBarry Smith   InsertMode  addv;
2421eb62cbbSBarry Smith   Scalar      *rvalues,*svalues;
2431eb62cbbSBarry Smith 
2441eb62cbbSBarry Smith   /* make sure all processors are either in INSERTMODE or ADDMODE */
24547794344SBarry Smith   MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
246dbb450caSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
247e3372554SBarry Smith     SETERRQ(1,0,"Some processors inserted others added");
2481eb62cbbSBarry Smith   }
24947794344SBarry Smith   mat->insertmode = addv; /* in case this processor had no cache */
2501eb62cbbSBarry Smith 
2511eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
2520452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
253cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
2540452661fSBarry Smith   owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
2551eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
2561eb62cbbSBarry Smith     idx = aij->stash.idx[i];
25717699dbbSLois Curfman McInnes     for ( j=0; j<size; j++ ) {
2581eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
2591eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
2608a729477SBarry Smith       }
2618a729477SBarry Smith     }
2628a729477SBarry Smith   }
26317699dbbSLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
2641eb62cbbSBarry Smith 
2651eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
2660452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
26717699dbbSLois Curfman McInnes   MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm);
26817699dbbSLois Curfman McInnes   nreceives = work[rank];
26917699dbbSLois Curfman McInnes   MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);
27017699dbbSLois Curfman McInnes   nmax = work[rank];
2710452661fSBarry Smith   PetscFree(work);
2721eb62cbbSBarry Smith 
2731eb62cbbSBarry Smith   /* post receives:
2741eb62cbbSBarry Smith        1) each message will consist of ordered pairs
2751eb62cbbSBarry Smith      (global index,value) we store the global index as a double
276d6dfbf8fSBarry Smith      to simplify the message passing.
2771eb62cbbSBarry Smith        2) since we don't know how long each individual message is we
2781eb62cbbSBarry Smith      allocate the largest needed buffer for each receive. Potentially
2791eb62cbbSBarry Smith      this is a lot of wasted space.
2801eb62cbbSBarry Smith 
2811eb62cbbSBarry Smith 
2821eb62cbbSBarry Smith        This could be done better.
2831eb62cbbSBarry Smith   */
2840452661fSBarry Smith   rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));
28578b31e54SBarry Smith   CHKPTRQ(rvalues);
2860452661fSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));
28778b31e54SBarry Smith   CHKPTRQ(recv_waits);
2881eb62cbbSBarry Smith   for ( i=0; i<nreceives; i++ ) {
289416022c9SBarry Smith     MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,
2901eb62cbbSBarry Smith               comm,recv_waits+i);
2911eb62cbbSBarry Smith   }
2921eb62cbbSBarry Smith 
2931eb62cbbSBarry Smith   /* do sends:
2941eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
2951eb62cbbSBarry Smith          the ith processor
2961eb62cbbSBarry Smith   */
2970452661fSBarry Smith   svalues = (Scalar *) PetscMalloc(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues);
2980452661fSBarry Smith   send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));
29978b31e54SBarry Smith   CHKPTRQ(send_waits);
3000452661fSBarry Smith   starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts);
3011eb62cbbSBarry Smith   starts[0] = 0;
30217699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3031eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
3041eb62cbbSBarry Smith     svalues[3*starts[owner[i]]]       = (Scalar)  aij->stash.idx[i];
3051eb62cbbSBarry Smith     svalues[3*starts[owner[i]]+1]     = (Scalar)  aij->stash.idy[i];
3061eb62cbbSBarry Smith     svalues[3*(starts[owner[i]]++)+2] =  aij->stash.array[i];
3071eb62cbbSBarry Smith   }
3080452661fSBarry Smith   PetscFree(owner);
3091eb62cbbSBarry Smith   starts[0] = 0;
31017699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3111eb62cbbSBarry Smith   count = 0;
31217699dbbSLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3131eb62cbbSBarry Smith     if (procs[i]) {
314416022c9SBarry Smith       MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,
3151eb62cbbSBarry Smith                 comm,send_waits+count++);
3161eb62cbbSBarry Smith     }
3171eb62cbbSBarry Smith   }
3180452661fSBarry Smith   PetscFree(starts); PetscFree(nprocs);
3191eb62cbbSBarry Smith 
3201eb62cbbSBarry Smith   /* Free cache space */
32155908cccSLois Curfman McInnes   PLogInfo(mat,"MatAssemblyBegin_MPIAIJ:Number of off-processor values %d\n",aij->stash.n);
32278b31e54SBarry Smith   ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr);
3231eb62cbbSBarry Smith 
3241eb62cbbSBarry Smith   aij->svalues    = svalues;    aij->rvalues    = rvalues;
3251eb62cbbSBarry Smith   aij->nsends     = nsends;     aij->nrecvs     = nreceives;
3261eb62cbbSBarry Smith   aij->send_waits = send_waits; aij->recv_waits = recv_waits;
3271eb62cbbSBarry Smith   aij->rmax       = nmax;
3281eb62cbbSBarry Smith 
3291eb62cbbSBarry Smith   return 0;
3301eb62cbbSBarry Smith }
33144a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat);
3321eb62cbbSBarry Smith 
3335615d1e5SSatish Balay #undef __FUNC__
3345615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIAIJ"
3358f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode)
3361eb62cbbSBarry Smith {
33744a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
3381eb62cbbSBarry Smith   MPI_Status  *send_status,recv_status;
339416022c9SBarry Smith   int         imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr;
340905e6a2fSBarry Smith   int         row,col,other_disassembled;
3411eb62cbbSBarry Smith   Scalar      *values,val;
34247794344SBarry Smith   InsertMode  addv = mat->insertmode;
3431eb62cbbSBarry Smith 
3441eb62cbbSBarry Smith   /*  wait on receives */
3451eb62cbbSBarry Smith   while (count) {
346d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);
3471eb62cbbSBarry Smith     /* unpack receives into our local space */
348d6dfbf8fSBarry Smith     values = aij->rvalues + 3*imdex*aij->rmax;
349c60448a5SBarry Smith     MPI_Get_count(&recv_status,MPIU_SCALAR,&n);
3501eb62cbbSBarry Smith     n = n/3;
3511eb62cbbSBarry Smith     for ( i=0; i<n; i++ ) {
352227d817aSBarry Smith       row = (int) PetscReal(values[3*i]) - aij->rstart;
353227d817aSBarry Smith       col = (int) PetscReal(values[3*i+1]);
3541eb62cbbSBarry Smith       val = values[3*i+2];
3551eb62cbbSBarry Smith       if (col >= aij->cstart && col < aij->cend) {
3561eb62cbbSBarry Smith         col -= aij->cstart;
3571eb62cbbSBarry Smith         MatSetValues(aij->A,1,&row,1,&col,&val,addv);
3581eb62cbbSBarry Smith       }
3591eb62cbbSBarry Smith       else {
36055a1b374SBarry Smith         if (mat->was_assembled || mat->assembled) {
361905e6a2fSBarry Smith           if (!aij->colmap) {
362905e6a2fSBarry Smith             ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
363905e6a2fSBarry Smith           }
364905e6a2fSBarry Smith           col = aij->colmap[col] - 1;
365ec8511deSBarry Smith           if (col < 0  && !((Mat_SeqAIJ*)(aij->A->data))->nonew) {
3662493cbb0SBarry Smith             ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
367227d817aSBarry Smith             col = (int) PetscReal(values[3*i+1]);
368d6dfbf8fSBarry Smith           }
3699e25ed09SBarry Smith         }
3701eb62cbbSBarry Smith         MatSetValues(aij->B,1,&row,1,&col,&val,addv);
3711eb62cbbSBarry Smith       }
3721eb62cbbSBarry Smith     }
3731eb62cbbSBarry Smith     count--;
3741eb62cbbSBarry Smith   }
3750452661fSBarry Smith   PetscFree(aij->recv_waits); PetscFree(aij->rvalues);
3761eb62cbbSBarry Smith 
3771eb62cbbSBarry Smith   /* wait on sends */
3781eb62cbbSBarry Smith   if (aij->nsends) {
3790a198c4cSBarry Smith     send_status = (MPI_Status *) PetscMalloc(aij->nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
3801eb62cbbSBarry Smith     MPI_Waitall(aij->nsends,aij->send_waits,send_status);
3810452661fSBarry Smith     PetscFree(send_status);
3821eb62cbbSBarry Smith   }
3830452661fSBarry Smith   PetscFree(aij->send_waits); PetscFree(aij->svalues);
3841eb62cbbSBarry Smith 
38578b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr);
38678b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr);
3871eb62cbbSBarry Smith 
3882493cbb0SBarry Smith   /* determine if any processor has disassembled, if so we must
3892493cbb0SBarry Smith      also disassemble ourselfs, in order that we may reassemble. */
390227d817aSBarry Smith   MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);
391227d817aSBarry Smith   if (mat->was_assembled && !other_disassembled) {
3922493cbb0SBarry Smith     ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
3932493cbb0SBarry Smith   }
3942493cbb0SBarry Smith 
3956d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
39678b31e54SBarry Smith     ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr);
3975e42470aSBarry Smith   }
39878b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr);
39978b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr);
4005e42470aSBarry Smith 
4017a0afa10SBarry Smith   if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;}
4028a729477SBarry Smith   return 0;
4038a729477SBarry Smith }
4048a729477SBarry Smith 
4055615d1e5SSatish Balay #undef __FUNC__
4065615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIAIJ"
4078f6be9afSLois Curfman McInnes int MatZeroEntries_MPIAIJ(Mat A)
4081eb62cbbSBarry Smith {
40944a69424SLois Curfman McInnes   Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data;
410dbd7a890SLois Curfman McInnes   int        ierr;
41178b31e54SBarry Smith   ierr = MatZeroEntries(l->A); CHKERRQ(ierr);
41278b31e54SBarry Smith   ierr = MatZeroEntries(l->B); CHKERRQ(ierr);
4131eb62cbbSBarry Smith   return 0;
4141eb62cbbSBarry Smith }
4151eb62cbbSBarry Smith 
4161eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the
4171eb62cbbSBarry Smith    matrix is square and the column and row owerships are identical.
4181eb62cbbSBarry Smith    This is a BUG. The only way to fix it seems to be to access
4191eb62cbbSBarry Smith    aij->A and aij->B directly and not through the MatZeroRows()
4201eb62cbbSBarry Smith    routine.
4211eb62cbbSBarry Smith */
4225615d1e5SSatish Balay #undef __FUNC__
4235615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ"
4248f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag)
4251eb62cbbSBarry Smith {
42644a69424SLois Curfman McInnes   Mat_MPIAIJ     *l = (Mat_MPIAIJ *) A->data;
42717699dbbSLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
4286abc6512SBarry Smith   int            *procs,*nprocs,j,found,idx,nsends,*work;
42917699dbbSLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
4305392566eSBarry Smith   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
431d6dfbf8fSBarry Smith   int            *lens,imdex,*lrows,*values;
432d6dfbf8fSBarry Smith   MPI_Comm       comm = A->comm;
4331eb62cbbSBarry Smith   MPI_Request    *send_waits,*recv_waits;
4341eb62cbbSBarry Smith   MPI_Status     recv_status,*send_status;
4351eb62cbbSBarry Smith   IS             istmp;
4361eb62cbbSBarry Smith 
43777c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
43878b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
4391eb62cbbSBarry Smith 
4401eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
4410452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
442cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
4430452661fSBarry Smith   owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
4441eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
4451eb62cbbSBarry Smith     idx = rows[i];
4461eb62cbbSBarry Smith     found = 0;
44717699dbbSLois Curfman McInnes     for ( j=0; j<size; j++ ) {
4481eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
4491eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
4501eb62cbbSBarry Smith       }
4511eb62cbbSBarry Smith     }
452e3372554SBarry Smith     if (!found) SETERRQ(1,0,"Index out of range");
4531eb62cbbSBarry Smith   }
45417699dbbSLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
4551eb62cbbSBarry Smith 
4561eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
4570452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
45817699dbbSLois Curfman McInnes   MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);
45917699dbbSLois Curfman McInnes   nrecvs = work[rank];
46017699dbbSLois Curfman McInnes   MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);
46117699dbbSLois Curfman McInnes   nmax = work[rank];
4620452661fSBarry Smith   PetscFree(work);
4631eb62cbbSBarry Smith 
4641eb62cbbSBarry Smith   /* post receives:   */
4650452661fSBarry Smith   rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */
46678b31e54SBarry Smith   CHKPTRQ(rvalues);
4670452661fSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));
46878b31e54SBarry Smith   CHKPTRQ(recv_waits);
4691eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
470416022c9SBarry Smith     MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
4711eb62cbbSBarry Smith   }
4721eb62cbbSBarry Smith 
4731eb62cbbSBarry Smith   /* do sends:
4741eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
4751eb62cbbSBarry Smith          the ith processor
4761eb62cbbSBarry Smith   */
4770452661fSBarry Smith   svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
4780452661fSBarry Smith   send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));
47978b31e54SBarry Smith   CHKPTRQ(send_waits);
4800452661fSBarry Smith   starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts);
4811eb62cbbSBarry Smith   starts[0] = 0;
48217699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
4831eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
4841eb62cbbSBarry Smith     svalues[starts[owner[i]]++] = rows[i];
4851eb62cbbSBarry Smith   }
4861eb62cbbSBarry Smith   ISRestoreIndices(is,&rows);
4871eb62cbbSBarry Smith 
4881eb62cbbSBarry Smith   starts[0] = 0;
48917699dbbSLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
4901eb62cbbSBarry Smith   count = 0;
49117699dbbSLois Curfman McInnes   for ( i=0; i<size; i++ ) {
4921eb62cbbSBarry Smith     if (procs[i]) {
493416022c9SBarry Smith       MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);
4941eb62cbbSBarry Smith     }
4951eb62cbbSBarry Smith   }
4960452661fSBarry Smith   PetscFree(starts);
4971eb62cbbSBarry Smith 
49817699dbbSLois Curfman McInnes   base = owners[rank];
4991eb62cbbSBarry Smith 
5001eb62cbbSBarry Smith   /*  wait on receives */
5010452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
5021eb62cbbSBarry Smith   source = lens + nrecvs;
5031eb62cbbSBarry Smith   count  = nrecvs; slen = 0;
5041eb62cbbSBarry Smith   while (count) {
505d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
5061eb62cbbSBarry Smith     /* unpack receives into our local space */
5071eb62cbbSBarry Smith     MPI_Get_count(&recv_status,MPI_INT,&n);
508d6dfbf8fSBarry Smith     source[imdex]  = recv_status.MPI_SOURCE;
509d6dfbf8fSBarry Smith     lens[imdex]  = n;
5101eb62cbbSBarry Smith     slen += n;
5111eb62cbbSBarry Smith     count--;
5121eb62cbbSBarry Smith   }
5130452661fSBarry Smith   PetscFree(recv_waits);
5141eb62cbbSBarry Smith 
5151eb62cbbSBarry Smith   /* move the data into the send scatter */
5160452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
5171eb62cbbSBarry Smith   count = 0;
5181eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
5191eb62cbbSBarry Smith     values = rvalues + i*nmax;
5201eb62cbbSBarry Smith     for ( j=0; j<lens[i]; j++ ) {
5211eb62cbbSBarry Smith       lrows[count++] = values[j] - base;
5221eb62cbbSBarry Smith     }
5231eb62cbbSBarry Smith   }
5240452661fSBarry Smith   PetscFree(rvalues); PetscFree(lens);
5250452661fSBarry Smith   PetscFree(owner); PetscFree(nprocs);
5261eb62cbbSBarry Smith 
5271eb62cbbSBarry Smith   /* actually zap the local rows */
528537820f0SBarry Smith   ierr = ISCreateGeneral(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
529464493b3SBarry Smith   PLogObjectParent(A,istmp);
5300452661fSBarry Smith   PetscFree(lrows);
53178b31e54SBarry Smith   ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
53278b31e54SBarry Smith   ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr);
53378b31e54SBarry Smith   ierr = ISDestroy(istmp); CHKERRQ(ierr);
5341eb62cbbSBarry Smith 
5351eb62cbbSBarry Smith   /* wait on sends */
5361eb62cbbSBarry Smith   if (nsends) {
5370452661fSBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));
53878b31e54SBarry Smith     CHKPTRQ(send_status);
5391eb62cbbSBarry Smith     MPI_Waitall(nsends,send_waits,send_status);
5400452661fSBarry Smith     PetscFree(send_status);
5411eb62cbbSBarry Smith   }
5420452661fSBarry Smith   PetscFree(send_waits); PetscFree(svalues);
5431eb62cbbSBarry Smith 
5441eb62cbbSBarry Smith   return 0;
5451eb62cbbSBarry Smith }
5461eb62cbbSBarry Smith 
5475615d1e5SSatish Balay #undef __FUNC__
5485615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ"
5498f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy)
5501eb62cbbSBarry Smith {
551416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
552fbd6ef76SBarry Smith   int        ierr,nt;
553416022c9SBarry Smith 
554a2ce50c7SBarry Smith   ierr = VecGetLocalSize(xx,&nt);  CHKERRQ(ierr);
555fbd6ef76SBarry Smith   if (nt != a->n) {
556f40265daSLois Curfman McInnes     SETERRQ(1,0,"Incompatible partition of A and xx");
557fbd6ef76SBarry Smith   }
55843a90d84SBarry Smith   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr);
55938597bd4SSatish Balay   ierr = (*a->A->ops.mult)(a->A,xx,yy); CHKERRQ(ierr);
56043a90d84SBarry Smith   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr);
56138597bd4SSatish Balay   ierr = (*a->B->ops.multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr);
5621eb62cbbSBarry Smith   return 0;
5631eb62cbbSBarry Smith }
5641eb62cbbSBarry Smith 
5655615d1e5SSatish Balay #undef __FUNC__
5665615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ"
5678f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
568da3a660dSBarry Smith {
569416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
570da3a660dSBarry Smith   int        ierr;
57143a90d84SBarry Smith   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
57238597bd4SSatish Balay   ierr = (*a->A->ops.multadd)(a->A,xx,yy,zz); CHKERRQ(ierr);
57343a90d84SBarry Smith   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
57438597bd4SSatish Balay   ierr = (*a->B->ops.multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr);
575da3a660dSBarry Smith   return 0;
576da3a660dSBarry Smith }
577da3a660dSBarry Smith 
5785615d1e5SSatish Balay #undef __FUNC__
5795615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ"
5808f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy)
581da3a660dSBarry Smith {
582416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
583da3a660dSBarry Smith   int        ierr;
584da3a660dSBarry Smith 
585da3a660dSBarry Smith   /* do nondiagonal part */
58638597bd4SSatish Balay   ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr);
587da3a660dSBarry Smith   /* send it on its way */
588537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
589da3a660dSBarry Smith   /* do local part */
59038597bd4SSatish Balay   ierr = (*a->A->ops.multtrans)(a->A,xx,yy); CHKERRQ(ierr);
591da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
592da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
593da3a660dSBarry Smith   /* but is not perhaps always true. */
594537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
595da3a660dSBarry Smith   return 0;
596da3a660dSBarry Smith }
597da3a660dSBarry Smith 
5985615d1e5SSatish Balay #undef __FUNC__
5995615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ"
6008f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
601da3a660dSBarry Smith {
602416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
603da3a660dSBarry Smith   int        ierr;
604da3a660dSBarry Smith 
605da3a660dSBarry Smith   /* do nondiagonal part */
60638597bd4SSatish Balay   ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr);
607da3a660dSBarry Smith   /* send it on its way */
608537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
609da3a660dSBarry Smith   /* do local part */
61038597bd4SSatish Balay   ierr = (*a->A->ops.multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr);
611da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
612da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
613da3a660dSBarry Smith   /* but is not perhaps always true. */
6140a198c4cSBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
615da3a660dSBarry Smith   return 0;
616da3a660dSBarry Smith }
617da3a660dSBarry Smith 
6181eb62cbbSBarry Smith /*
6191eb62cbbSBarry Smith   This only works correctly for square matrices where the subblock A->A is the
6201eb62cbbSBarry Smith    diagonal block
6211eb62cbbSBarry Smith */
6225615d1e5SSatish Balay #undef __FUNC__
6235615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ"
6248f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v)
6251eb62cbbSBarry Smith {
626416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
62744cd7ae7SLois Curfman McInnes   if (a->M != a->N)
628e3372554SBarry Smith     SETERRQ(1,0,"Supports only square matrix where A->A is diag block");
6295baf8537SBarry Smith   if (a->rstart != a->cstart || a->rend != a->cend) {
630e3372554SBarry Smith     SETERRQ(1,0,"row partition must equal col partition");  }
631416022c9SBarry Smith   return MatGetDiagonal(a->A,v);
6321eb62cbbSBarry Smith }
6331eb62cbbSBarry Smith 
6345615d1e5SSatish Balay #undef __FUNC__
6355615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ"
6368f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A)
637052efed2SBarry Smith {
638052efed2SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
639052efed2SBarry Smith   int        ierr;
640052efed2SBarry Smith   ierr = MatScale(aa,a->A); CHKERRQ(ierr);
641052efed2SBarry Smith   ierr = MatScale(aa,a->B); CHKERRQ(ierr);
642052efed2SBarry Smith   return 0;
643052efed2SBarry Smith }
644052efed2SBarry Smith 
6455615d1e5SSatish Balay #undef __FUNC__
6465eea60f9SBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ" /* ADIC Ignore */
6478f6be9afSLois Curfman McInnes int MatDestroy_MPIAIJ(PetscObject obj)
6481eb62cbbSBarry Smith {
6491eb62cbbSBarry Smith   Mat        mat = (Mat) obj;
65044a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
6511eb62cbbSBarry Smith   int        ierr;
652a5a9c739SBarry Smith #if defined(PETSC_LOG)
6536e35fa57SLois Curfman McInnes   PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N);
654a5a9c739SBarry Smith #endif
6550452661fSBarry Smith   PetscFree(aij->rowners);
65678b31e54SBarry Smith   ierr = MatDestroy(aij->A); CHKERRQ(ierr);
65778b31e54SBarry Smith   ierr = MatDestroy(aij->B); CHKERRQ(ierr);
6580452661fSBarry Smith   if (aij->colmap) PetscFree(aij->colmap);
6590452661fSBarry Smith   if (aij->garray) PetscFree(aij->garray);
6601eb62cbbSBarry Smith   if (aij->lvec)   VecDestroy(aij->lvec);
661a56f8943SBarry Smith   if (aij->Mvctx)  VecScatterDestroy(aij->Mvctx);
6627a0afa10SBarry Smith   if (aij->rowvalues) PetscFree(aij->rowvalues);
6630452661fSBarry Smith   PetscFree(aij);
66490f02eecSBarry Smith   if (mat->mapping) {
66590f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
66690f02eecSBarry Smith   }
667a5a9c739SBarry Smith   PLogObjectDestroy(mat);
6680452661fSBarry Smith   PetscHeaderDestroy(mat);
6691eb62cbbSBarry Smith   return 0;
6701eb62cbbSBarry Smith }
671ee50ffe9SBarry Smith 
6725615d1e5SSatish Balay #undef __FUNC__
6735eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary" /* ADIC Ignore */
6748f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer)
6751eb62cbbSBarry Smith {
676416022c9SBarry Smith   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
677416022c9SBarry Smith   int         ierr;
678416022c9SBarry Smith 
67917699dbbSLois Curfman McInnes   if (aij->size == 1) {
680416022c9SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
681416022c9SBarry Smith   }
682e3372554SBarry Smith   else SETERRQ(1,0,"Only uniprocessor output supported");
683416022c9SBarry Smith   return 0;
684416022c9SBarry Smith }
685416022c9SBarry Smith 
6865615d1e5SSatish Balay #undef __FUNC__
6875eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworMatlab" /* ADIC Ignore */
6888f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer)
689416022c9SBarry Smith {
69044a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
691dbb450caSBarry Smith   Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data;
692a56f8943SBarry Smith   int         ierr, format,shift = C->indexshift,rank;
693d636dbe3SBarry Smith   FILE        *fd;
69419bcc07fSBarry Smith   ViewerType  vtype;
695416022c9SBarry Smith 
69619bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
69719bcc07fSBarry Smith   if (vtype  == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) {
69890ace30eSBarry Smith     ierr = ViewerGetFormat(viewer,&format);
6990a198c4cSBarry Smith     if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
7004e220ebcSLois Curfman McInnes       MatInfo info;
7014e220ebcSLois Curfman McInnes       int     flg;
702a56f8943SBarry Smith       MPI_Comm_rank(mat->comm,&rank);
70390ace30eSBarry Smith       ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
7044e220ebcSLois Curfman McInnes       ierr = MatGetInfo(mat,MAT_LOCAL,&info);
70595e01e2fSLois Curfman McInnes       ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr);
70677c4ece6SBarry Smith       PetscSequentialPhaseBegin(mat->comm,1);
70795e01e2fSLois Curfman McInnes       if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n",
7084e220ebcSLois Curfman McInnes          rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
70995e01e2fSLois Curfman McInnes       else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n",
7104e220ebcSLois Curfman McInnes          rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
7114e220ebcSLois Curfman McInnes       ierr = MatGetInfo(aij->A,MAT_LOCAL,&info);
7124e220ebcSLois Curfman McInnes       fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used);
7134e220ebcSLois Curfman McInnes       ierr = MatGetInfo(aij->B,MAT_LOCAL,&info);
7144e220ebcSLois Curfman McInnes       fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used);
715a56f8943SBarry Smith       fflush(fd);
71677c4ece6SBarry Smith       PetscSequentialPhaseEnd(mat->comm,1);
717a40aa06bSLois Curfman McInnes       ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr);
718416022c9SBarry Smith       return 0;
719416022c9SBarry Smith     }
7200a198c4cSBarry Smith     else if (format == VIEWER_FORMAT_ASCII_INFO) {
72108480c60SBarry Smith       return 0;
72208480c60SBarry Smith     }
723416022c9SBarry Smith   }
724416022c9SBarry Smith 
72519bcc07fSBarry Smith   if (vtype == DRAW_VIEWER) {
72619bcc07fSBarry Smith     Draw       draw;
72719bcc07fSBarry Smith     PetscTruth isnull;
72819bcc07fSBarry Smith     ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
72919bcc07fSBarry Smith     ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
73019bcc07fSBarry Smith   }
73119bcc07fSBarry Smith 
73219bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
73390ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
73477c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
735d13ab20cSBarry Smith     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n",
73617699dbbSLois Curfman McInnes            aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart,
7371eb62cbbSBarry Smith            aij->cend);
73878b31e54SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
73978b31e54SBarry Smith     ierr = MatView(aij->B,viewer); CHKERRQ(ierr);
740d13ab20cSBarry Smith     fflush(fd);
74177c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
742d13ab20cSBarry Smith   }
743416022c9SBarry Smith   else {
744a56f8943SBarry Smith     int size = aij->size;
745a56f8943SBarry Smith     rank = aij->rank;
74617699dbbSLois Curfman McInnes     if (size == 1) {
74778b31e54SBarry Smith       ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
74895373324SBarry Smith     }
74995373324SBarry Smith     else {
75095373324SBarry Smith       /* assemble the entire matrix onto first processor. */
75195373324SBarry Smith       Mat         A;
752ec8511deSBarry Smith       Mat_SeqAIJ *Aloc;
7532eb8c8abSBarry Smith       int         M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct;
75495373324SBarry Smith       Scalar      *a;
7552ee70a88SLois Curfman McInnes 
75617699dbbSLois Curfman McInnes       if (!rank) {
757b4fd4287SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);
758c451ab25SLois Curfman McInnes                CHKERRQ(ierr);
75995373324SBarry Smith       }
76095373324SBarry Smith       else {
761b4fd4287SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);
762c451ab25SLois Curfman McInnes                CHKERRQ(ierr);
76395373324SBarry Smith       }
764464493b3SBarry Smith       PLogObjectParent(mat,A);
765416022c9SBarry Smith 
76695373324SBarry Smith       /* copy over the A part */
767ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->A->data;
7682ee70a88SLois Curfman McInnes       m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
76995373324SBarry Smith       row = aij->rstart;
770dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;}
77195373324SBarry Smith       for ( i=0; i<m; i++ ) {
772416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr);
77395373324SBarry Smith         row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
77495373324SBarry Smith       }
7752ee70a88SLois Curfman McInnes       aj = Aloc->j;
776dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;}
77795373324SBarry Smith 
77895373324SBarry Smith       /* copy over the B part */
779ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->B->data;
7802ee70a88SLois Curfman McInnes       m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
78195373324SBarry Smith       row = aij->rstart;
7820452661fSBarry Smith       ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
783dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];}
78495373324SBarry Smith       for ( i=0; i<m; i++ ) {
785416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr);
78695373324SBarry Smith         row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
78795373324SBarry Smith       }
7880452661fSBarry Smith       PetscFree(ct);
7896d4a8577SBarry Smith       ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
7906d4a8577SBarry Smith       ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
79117699dbbSLois Curfman McInnes       if (!rank) {
79278b31e54SBarry Smith         ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr);
79395373324SBarry Smith       }
79478b31e54SBarry Smith       ierr = MatDestroy(A); CHKERRQ(ierr);
79595373324SBarry Smith     }
79695373324SBarry Smith   }
7971eb62cbbSBarry Smith   return 0;
7981eb62cbbSBarry Smith }
7991eb62cbbSBarry Smith 
8005615d1e5SSatish Balay #undef __FUNC__
8015eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ" /* ADIC Ignore */
8028f6be9afSLois Curfman McInnes int MatView_MPIAIJ(PetscObject obj,Viewer viewer)
803416022c9SBarry Smith {
804416022c9SBarry Smith   Mat         mat = (Mat) obj;
805416022c9SBarry Smith   int         ierr;
80619bcc07fSBarry Smith   ViewerType  vtype;
807416022c9SBarry Smith 
80819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
80919bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER ||
81019bcc07fSBarry Smith       vtype == DRAW_VIEWER       || vtype == MATLAB_VIEWER) {
811d7e8b826SBarry Smith     ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr);
812416022c9SBarry Smith   }
81319bcc07fSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
814416022c9SBarry Smith     return MatView_MPIAIJ_Binary(mat,viewer);
815416022c9SBarry Smith   }
816416022c9SBarry Smith   return 0;
817416022c9SBarry Smith }
818416022c9SBarry Smith 
8191eb62cbbSBarry Smith /*
8201eb62cbbSBarry Smith     This has to provide several versions.
8211eb62cbbSBarry Smith 
8221eb62cbbSBarry Smith      2) a) use only local smoothing updating outer values only once.
8231eb62cbbSBarry Smith         b) local smoothing updating outer values each inner iteration
824d6dfbf8fSBarry Smith      3) color updating out values betwen colors.
8251eb62cbbSBarry Smith */
8265615d1e5SSatish Balay #undef __FUNC__
8275615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ"
8288f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag,
829dbb450caSBarry Smith                            double fshift,int its,Vec xx)
8308a729477SBarry Smith {
83144a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
832d6dfbf8fSBarry Smith   Mat        AA = mat->A, BB = mat->B;
833ec8511deSBarry Smith   Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data;
834c16cb8f2SBarry Smith   Scalar     *b,*x,*xs,*ls,d,*v,sum;
8356abc6512SBarry Smith   int        ierr,*idx, *diag;
836416022c9SBarry Smith   int        n = mat->n, m = mat->m, i,shift = A->indexshift;
8378a729477SBarry Smith 
838d6dfbf8fSBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls);
839dbb450caSBarry Smith   xs = x + shift; /* shift by one for index start of 1 */
840dbb450caSBarry Smith   ls = ls + shift;
841ec8511deSBarry Smith   if (!A->diag) {if ((ierr = MatMarkDiag_SeqAIJ(AA))) return ierr;}
842d6dfbf8fSBarry Smith   diag = A->diag;
843c16cb8f2SBarry Smith   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
844da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
84538597bd4SSatish Balay       return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);
846da3a660dSBarry Smith     }
84743a90d84SBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
84878b31e54SBarry Smith     CHKERRQ(ierr);
84943a90d84SBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
85078b31e54SBarry Smith     CHKERRQ(ierr);
851d6dfbf8fSBarry Smith     while (its--) {
852d6dfbf8fSBarry Smith       /* go down through the rows */
853d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
854d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
855dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
856dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
857d6dfbf8fSBarry Smith         sum  = b[i];
858d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
859dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
860d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
861dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
862dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
863d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
86455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
865d6dfbf8fSBarry Smith       }
866d6dfbf8fSBarry Smith       /* come up through the rows */
867d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
868d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
869dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
870dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
871d6dfbf8fSBarry Smith         sum  = b[i];
872d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
873dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
874d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
875dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
876dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
877d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
87855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
879d6dfbf8fSBarry Smith       }
880d6dfbf8fSBarry Smith     }
881d6dfbf8fSBarry Smith   }
882d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_FORWARD_SWEEP){
883da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
88438597bd4SSatish Balay       return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);
885da3a660dSBarry Smith     }
88643a90d84SBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
88778b31e54SBarry Smith     CHKERRQ(ierr);
88843a90d84SBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
88978b31e54SBarry Smith     CHKERRQ(ierr);
890d6dfbf8fSBarry Smith     while (its--) {
891d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
892d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
893dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
894dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
895d6dfbf8fSBarry Smith         sum  = b[i];
896d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
897dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
898d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
899dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
900dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
901d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
90255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
903d6dfbf8fSBarry Smith       }
904d6dfbf8fSBarry Smith     }
905d6dfbf8fSBarry Smith   }
906d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_BACKWARD_SWEEP){
907da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
90838597bd4SSatish Balay       return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);
909da3a660dSBarry Smith     }
91043a90d84SBarry Smith     ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,
91178b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
91243a90d84SBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,
91378b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
914d6dfbf8fSBarry Smith     while (its--) {
915d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
916d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
917dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
918dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
919d6dfbf8fSBarry Smith         sum  = b[i];
920d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
921dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
922d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
923dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
924dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
925d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
92655a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
927d6dfbf8fSBarry Smith       }
928d6dfbf8fSBarry Smith     }
929d6dfbf8fSBarry Smith   }
930c16cb8f2SBarry Smith   else {
931e3372554SBarry Smith     SETERRQ(1,0,"Option not supported");
932c16cb8f2SBarry Smith   }
9338a729477SBarry Smith   return 0;
9348a729477SBarry Smith }
935a66be287SLois Curfman McInnes 
9365615d1e5SSatish Balay #undef __FUNC__
9375eea60f9SBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ" /* ADIC Ignore */
9388f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info)
939a66be287SLois Curfman McInnes {
940a66be287SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
941a66be287SLois Curfman McInnes   Mat        A = mat->A, B = mat->B;
9424e220ebcSLois Curfman McInnes   int        ierr;
9434e220ebcSLois Curfman McInnes   double     isend[5], irecv[5];
944a66be287SLois Curfman McInnes 
9454e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
9464e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
9474e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
9484e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
9494e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
9504e220ebcSLois Curfman McInnes   ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr);
9514e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
9524e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
9534e220ebcSLois Curfman McInnes   ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr);
9544e220ebcSLois Curfman McInnes   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
9554e220ebcSLois Curfman McInnes   isend[3] += info->memory;  isend[4] += info->mallocs;
956a66be287SLois Curfman McInnes   if (flag == MAT_LOCAL) {
9574e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
9584e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
9594e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
9604e220ebcSLois Curfman McInnes     info->memory       = isend[3];
9614e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
962a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
9634e220ebcSLois Curfman McInnes     MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm);
9644e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
9654e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
9664e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
9674e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
9684e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
969a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
9704e220ebcSLois Curfman McInnes     MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm);
9714e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
9724e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
9734e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
9744e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
9754e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
976a66be287SLois Curfman McInnes   }
9774e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
9784e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
9794e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
9804e220ebcSLois Curfman McInnes 
981a66be287SLois Curfman McInnes   return 0;
982a66be287SLois Curfman McInnes }
983a66be287SLois Curfman McInnes 
984299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*);
985299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*);
986299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double);
987299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *);
988299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec);
989299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec);
990299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec);
991299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec);
992299609e3SLois Curfman McInnes 
9935615d1e5SSatish Balay #undef __FUNC__
9945eea60f9SBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ" /* ADIC Ignore */
9958f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op)
996c74985f6SBarry Smith {
997c0bbcb79SLois Curfman McInnes   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
998c74985f6SBarry Smith 
9996d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
10006d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
10016da5968aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED ||
1002c2653b3dSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
1003c2653b3dSLois Curfman McInnes       op == MAT_NEW_NONZERO_LOCATION_ERROR) {
1004b1fbbac0SLois Curfman McInnes         MatSetOption(a->A,op);
1005b1fbbac0SLois Curfman McInnes         MatSetOption(a->B,op);
1006b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
1007aeafbbfcSLois Curfman McInnes     a->roworiented = 1;
1008c0bbcb79SLois Curfman McInnes     MatSetOption(a->A,op);
1009c0bbcb79SLois Curfman McInnes     MatSetOption(a->B,op);
1010b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
10116da5968aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
10126d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
10136d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
10146d4a8577SBarry Smith              op == MAT_YES_NEW_DIAGONALS)
101594a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_MPIAIJ:Option ignored\n");
10166d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED) {
10174b0e389bSBarry Smith     a->roworiented = 0;
10184b0e389bSBarry Smith     MatSetOption(a->A,op);
10194b0e389bSBarry Smith     MatSetOption(a->B,op);
102090f02eecSBarry Smith   } else if (op == MAT_IGNORE_OFF_PROCESSOR_ENTRIES) {
102190f02eecSBarry Smith     a->donotstash = 1;
102290f02eecSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS)
1023e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");}
1024c0bbcb79SLois Curfman McInnes   else
1025e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
1026c74985f6SBarry Smith   return 0;
1027c74985f6SBarry Smith }
1028c74985f6SBarry Smith 
10295615d1e5SSatish Balay #undef __FUNC__
10305eea60f9SBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ" /* ADIC Ignore */
10318f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n)
1032c74985f6SBarry Smith {
103344a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1034c74985f6SBarry Smith   *m = mat->M; *n = mat->N;
1035c74985f6SBarry Smith   return 0;
1036c74985f6SBarry Smith }
1037c74985f6SBarry Smith 
10385615d1e5SSatish Balay #undef __FUNC__
10395eea60f9SBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ" /* ADIC Ignore */
10408f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n)
1041c74985f6SBarry Smith {
104244a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1043b7c46309SBarry Smith   *m = mat->m; *n = mat->N;
1044c74985f6SBarry Smith   return 0;
1045c74985f6SBarry Smith }
1046c74985f6SBarry Smith 
10475615d1e5SSatish Balay #undef __FUNC__
10485eea60f9SBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ" /* ADIC Ignore */
10498f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n)
1050c74985f6SBarry Smith {
105144a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1052c74985f6SBarry Smith   *m = mat->rstart; *n = mat->rend;
1053c74985f6SBarry Smith   return 0;
1054c74985f6SBarry Smith }
1055c74985f6SBarry Smith 
10566d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**);
10576d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**);
10586d84be18SBarry Smith 
10595615d1e5SSatish Balay #undef __FUNC__
10605615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ"
10616d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v)
106239e00950SLois Curfman McInnes {
1063154123eaSLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
106470f0671dSBarry Smith   Scalar     *vworkA, *vworkB, **pvA, **pvB,*v_p;
1065154123eaSLois Curfman McInnes   int        i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart;
1066154123eaSLois Curfman McInnes   int        nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend;
106770f0671dSBarry Smith   int        *cmap, *idx_p;
106839e00950SLois Curfman McInnes 
1069e3372554SBarry Smith   if (mat->getrowactive == PETSC_TRUE) SETERRQ(1,0,"Already active");
10707a0afa10SBarry Smith   mat->getrowactive = PETSC_TRUE;
10717a0afa10SBarry Smith 
107270f0671dSBarry Smith   if (!mat->rowvalues && (idx || v)) {
10737a0afa10SBarry Smith     /*
10747a0afa10SBarry Smith         allocate enough space to hold information from the longest row.
10757a0afa10SBarry Smith     */
10767a0afa10SBarry Smith     Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data;
1077c16cb8f2SBarry Smith     int     max = 1,m = mat->m,tmp;
1078c16cb8f2SBarry Smith     for ( i=0; i<m; i++ ) {
10797a0afa10SBarry Smith       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i];
10807a0afa10SBarry Smith       if (max < tmp) { max = tmp; }
10817a0afa10SBarry Smith     }
10827a0afa10SBarry Smith     mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar)));
10837a0afa10SBarry Smith     CHKPTRQ(mat->rowvalues);
10847a0afa10SBarry Smith     mat->rowindices = (int *) (mat->rowvalues + max);
10857a0afa10SBarry Smith   }
10867a0afa10SBarry Smith 
1087e3372554SBarry Smith   if (row < rstart || row >= rend) SETERRQ(1,0,"Only local rows")
1088abc0e9e4SLois Curfman McInnes   lrow = row - rstart;
108939e00950SLois Curfman McInnes 
1090154123eaSLois Curfman McInnes   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1091154123eaSLois Curfman McInnes   if (!v)   {pvA = 0; pvB = 0;}
1092154123eaSLois Curfman McInnes   if (!idx) {pcA = 0; if (!v) pcB = 0;}
109338597bd4SSatish Balay   ierr = (*mat->A->ops.getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
109438597bd4SSatish Balay   ierr = (*mat->B->ops.getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
1095154123eaSLois Curfman McInnes   nztot = nzA + nzB;
1096154123eaSLois Curfman McInnes 
109770f0671dSBarry Smith   cmap  = mat->garray;
1098154123eaSLois Curfman McInnes   if (v  || idx) {
1099154123eaSLois Curfman McInnes     if (nztot) {
1100154123eaSLois Curfman McInnes       /* Sort by increasing column numbers, assuming A and B already sorted */
110170f0671dSBarry Smith       int imark = -1;
1102154123eaSLois Curfman McInnes       if (v) {
110370f0671dSBarry Smith         *v = v_p = mat->rowvalues;
110439e00950SLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
110570f0671dSBarry Smith           if (cmap[cworkB[i]] < cstart)   v_p[i] = vworkB[i];
1106154123eaSLois Curfman McInnes           else break;
1107154123eaSLois Curfman McInnes         }
1108154123eaSLois Curfman McInnes         imark = i;
110970f0671dSBarry Smith         for ( i=0; i<nzA; i++ )     v_p[imark+i] = vworkA[i];
111070f0671dSBarry Smith         for ( i=imark; i<nzB; i++ ) v_p[nzA+i]   = vworkB[i];
1111154123eaSLois Curfman McInnes       }
1112154123eaSLois Curfman McInnes       if (idx) {
111370f0671dSBarry Smith         *idx = idx_p = mat->rowindices;
111470f0671dSBarry Smith         if (imark > -1) {
111570f0671dSBarry Smith           for ( i=0; i<imark; i++ ) {
111670f0671dSBarry Smith             idx_p[i] = cmap[cworkB[i]];
111770f0671dSBarry Smith           }
111870f0671dSBarry Smith         } else {
1119154123eaSLois Curfman McInnes           for ( i=0; i<nzB; i++ ) {
112070f0671dSBarry Smith             if (cmap[cworkB[i]] < cstart)   idx_p[i] = cmap[cworkB[i]];
1121154123eaSLois Curfman McInnes             else break;
1122154123eaSLois Curfman McInnes           }
1123154123eaSLois Curfman McInnes           imark = i;
112470f0671dSBarry Smith         }
112570f0671dSBarry Smith         for ( i=0; i<nzA; i++ )     idx_p[imark+i] = cstart + cworkA[i];
112670f0671dSBarry Smith         for ( i=imark; i<nzB; i++ ) idx_p[nzA+i]   = cmap[cworkB[i]];
112739e00950SLois Curfman McInnes       }
112839e00950SLois Curfman McInnes     }
11291ca473b0SSatish Balay     else {
11301ca473b0SSatish Balay       if (idx) *idx = 0;
11311ca473b0SSatish Balay       if (v)   *v   = 0;
11321ca473b0SSatish Balay     }
1133154123eaSLois Curfman McInnes   }
113439e00950SLois Curfman McInnes   *nz = nztot;
113538597bd4SSatish Balay   ierr = (*mat->A->ops.restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
113638597bd4SSatish Balay   ierr = (*mat->B->ops.restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
113739e00950SLois Curfman McInnes   return 0;
113839e00950SLois Curfman McInnes }
113939e00950SLois Curfman McInnes 
11405615d1e5SSatish Balay #undef __FUNC__
11415eea60f9SBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ" /* ADIC Ignore */
11426d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v)
114339e00950SLois Curfman McInnes {
11447a0afa10SBarry Smith   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
11457a0afa10SBarry Smith   if (aij->getrowactive == PETSC_FALSE) {
1146e3372554SBarry Smith     SETERRQ(1,0,"MatGetRow not called");
11477a0afa10SBarry Smith   }
11487a0afa10SBarry Smith   aij->getrowactive = PETSC_FALSE;
114939e00950SLois Curfman McInnes   return 0;
115039e00950SLois Curfman McInnes }
115139e00950SLois Curfman McInnes 
11525615d1e5SSatish Balay #undef __FUNC__
11535615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ"
11548f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm)
1155855ac2c5SLois Curfman McInnes {
1156855ac2c5SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1157ec8511deSBarry Smith   Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data;
1158416022c9SBarry Smith   int        ierr, i, j, cstart = aij->cstart,shift = amat->indexshift;
1159416022c9SBarry Smith   double     sum = 0.0;
116004ca555eSLois Curfman McInnes   Scalar     *v;
116104ca555eSLois Curfman McInnes 
116217699dbbSLois Curfman McInnes   if (aij->size == 1) {
116314183eadSLois Curfman McInnes     ierr =  MatNorm(aij->A,type,norm); CHKERRQ(ierr);
116437fa93a5SLois Curfman McInnes   } else {
116504ca555eSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
116604ca555eSLois Curfman McInnes       v = amat->a;
116704ca555eSLois Curfman McInnes       for (i=0; i<amat->nz; i++ ) {
116804ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
116904ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
117004ca555eSLois Curfman McInnes #else
117104ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
117204ca555eSLois Curfman McInnes #endif
117304ca555eSLois Curfman McInnes       }
117404ca555eSLois Curfman McInnes       v = bmat->a;
117504ca555eSLois Curfman McInnes       for (i=0; i<bmat->nz; i++ ) {
117604ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
117704ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
117804ca555eSLois Curfman McInnes #else
117904ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
118004ca555eSLois Curfman McInnes #endif
118104ca555eSLois Curfman McInnes       }
11826d84be18SBarry Smith       MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);
118304ca555eSLois Curfman McInnes       *norm = sqrt(*norm);
118404ca555eSLois Curfman McInnes     }
118504ca555eSLois Curfman McInnes     else if (type == NORM_1) { /* max column norm */
118604ca555eSLois Curfman McInnes       double *tmp, *tmp2;
118704ca555eSLois Curfman McInnes       int    *jj, *garray = aij->garray;
11880452661fSBarry Smith       tmp  = (double *) PetscMalloc( aij->N*sizeof(double) ); CHKPTRQ(tmp);
11890452661fSBarry Smith       tmp2 = (double *) PetscMalloc( aij->N*sizeof(double) ); CHKPTRQ(tmp2);
1190cddf8d76SBarry Smith       PetscMemzero(tmp,aij->N*sizeof(double));
119104ca555eSLois Curfman McInnes       *norm = 0.0;
119204ca555eSLois Curfman McInnes       v = amat->a; jj = amat->j;
119304ca555eSLois Curfman McInnes       for ( j=0; j<amat->nz; j++ ) {
1194579c6b6fSBarry Smith         tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v);  v++;
119504ca555eSLois Curfman McInnes       }
119604ca555eSLois Curfman McInnes       v = bmat->a; jj = bmat->j;
119704ca555eSLois Curfman McInnes       for ( j=0; j<bmat->nz; j++ ) {
1198579c6b6fSBarry Smith         tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++;
119904ca555eSLois Curfman McInnes       }
12006d84be18SBarry Smith       MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);
120104ca555eSLois Curfman McInnes       for ( j=0; j<aij->N; j++ ) {
120204ca555eSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
120304ca555eSLois Curfman McInnes       }
12040452661fSBarry Smith       PetscFree(tmp); PetscFree(tmp2);
120504ca555eSLois Curfman McInnes     }
120604ca555eSLois Curfman McInnes     else if (type == NORM_INFINITY) { /* max row norm */
1207515d9167SLois Curfman McInnes       double ntemp = 0.0;
120804ca555eSLois Curfman McInnes       for ( j=0; j<amat->m; j++ ) {
1209dbb450caSBarry Smith         v = amat->a + amat->i[j] + shift;
121004ca555eSLois Curfman McInnes         sum = 0.0;
121104ca555eSLois Curfman McInnes         for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) {
1212cddf8d76SBarry Smith           sum += PetscAbsScalar(*v); v++;
121304ca555eSLois Curfman McInnes         }
1214dbb450caSBarry Smith         v = bmat->a + bmat->i[j] + shift;
121504ca555eSLois Curfman McInnes         for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) {
1216cddf8d76SBarry Smith           sum += PetscAbsScalar(*v); v++;
121704ca555eSLois Curfman McInnes         }
1218515d9167SLois Curfman McInnes         if (sum > ntemp) ntemp = sum;
121904ca555eSLois Curfman McInnes       }
12206d84be18SBarry Smith       MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);
122104ca555eSLois Curfman McInnes     }
122204ca555eSLois Curfman McInnes     else {
1223e3372554SBarry Smith       SETERRQ(1,0,"No support for two norm");
122404ca555eSLois Curfman McInnes     }
122537fa93a5SLois Curfman McInnes   }
1226855ac2c5SLois Curfman McInnes   return 0;
1227855ac2c5SLois Curfman McInnes }
1228855ac2c5SLois Curfman McInnes 
12295615d1e5SSatish Balay #undef __FUNC__
12305615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ"
12318f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout)
1232b7c46309SBarry Smith {
1233b7c46309SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1234dbb450caSBarry Smith   Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data;
1235416022c9SBarry Smith   int        ierr,shift = Aloc->indexshift;
1236416022c9SBarry Smith   Mat        B;
1237b7c46309SBarry Smith   int        M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct;
1238b7c46309SBarry Smith   Scalar     *array;
1239b7c46309SBarry Smith 
12403638b69dSLois Curfman McInnes   if (matout == PETSC_NULL && M != N)
1241e3372554SBarry Smith     SETERRQ(1,0,"Square matrix only for in-place");
1242b4fd4287SBarry Smith   ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,PETSC_NULL,0,
1243b4fd4287SBarry Smith          PETSC_NULL,&B); CHKERRQ(ierr);
1244b7c46309SBarry Smith 
1245b7c46309SBarry Smith   /* copy over the A part */
1246ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->A->data;
1247b7c46309SBarry Smith   m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1248b7c46309SBarry Smith   row = a->rstart;
1249dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;}
1250b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1251416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1252b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
1253b7c46309SBarry Smith   }
1254b7c46309SBarry Smith   aj = Aloc->j;
12554af08d9eSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;}
1256b7c46309SBarry Smith 
1257b7c46309SBarry Smith   /* copy over the B part */
1258ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->B->data;
1259b7c46309SBarry Smith   m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1260b7c46309SBarry Smith   row = a->rstart;
12610452661fSBarry Smith   ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols);
1262dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];}
1263b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1264416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1265b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
1266b7c46309SBarry Smith   }
12670452661fSBarry Smith   PetscFree(ct);
12686d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12696d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12703638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
12710de55854SLois Curfman McInnes     *matout = B;
12720de55854SLois Curfman McInnes   } else {
12730de55854SLois Curfman McInnes     /* This isn't really an in-place transpose .... but free data structures from a */
12740452661fSBarry Smith     PetscFree(a->rowners);
12750de55854SLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
12760de55854SLois Curfman McInnes     ierr = MatDestroy(a->B); CHKERRQ(ierr);
12770452661fSBarry Smith     if (a->colmap) PetscFree(a->colmap);
12780452661fSBarry Smith     if (a->garray) PetscFree(a->garray);
12790de55854SLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
1280a56f8943SBarry Smith     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
12810452661fSBarry Smith     PetscFree(a);
1282416022c9SBarry Smith     PetscMemcpy(A,B,sizeof(struct _Mat));
12830452661fSBarry Smith     PetscHeaderDestroy(B);
12840de55854SLois Curfman McInnes   }
1285b7c46309SBarry Smith   return 0;
1286b7c46309SBarry Smith }
1287b7c46309SBarry Smith 
12885615d1e5SSatish Balay #undef __FUNC__
12895615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ"
12904b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr)
1291a008b906SSatish Balay {
12924b967eb1SSatish Balay   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
12934b967eb1SSatish Balay   Mat a = aij->A, b = aij->B;
1294a008b906SSatish Balay   int ierr,s1,s2,s3;
1295a008b906SSatish Balay 
12964b967eb1SSatish Balay   ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr);
12974b967eb1SSatish Balay   if (rr) {
12984b967eb1SSatish Balay     s3 = aij->n;
12994b967eb1SSatish Balay     VecGetLocalSize_Fast(rr,s1);
1300e3372554SBarry Smith     if (s1!=s3) SETERRQ(1,0,"right vector non-conforming local size");
13014b967eb1SSatish Balay     /* Overlap communication with computation. */
130243a90d84SBarry Smith     ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr);
1303a008b906SSatish Balay   }
13044b967eb1SSatish Balay   if (ll) {
13054b967eb1SSatish Balay     VecGetLocalSize_Fast(ll,s1);
1306e3372554SBarry Smith     if (s1!=s2) SETERRQ(1,0,"left vector non-conforming local size");
1307c351683dSSatish Balay     ierr = (*b->ops.diagonalscale)(b,ll,0); CHKERRQ(ierr);
13084b967eb1SSatish Balay   }
13094b967eb1SSatish Balay   /* scale  the diagonal block */
1310c351683dSSatish Balay   ierr = (*a->ops.diagonalscale)(a,ll,rr); CHKERRQ(ierr);
13114b967eb1SSatish Balay 
13124b967eb1SSatish Balay   if (rr) {
13134b967eb1SSatish Balay     /* Do a scatter end and then right scale the off-diagonal block */
131443a90d84SBarry Smith     ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr);
1315c351683dSSatish Balay     ierr = (*b->ops.diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr);
13164b967eb1SSatish Balay   }
13174b967eb1SSatish Balay 
1318a008b906SSatish Balay   return 0;
1319a008b906SSatish Balay }
1320a008b906SSatish Balay 
1321a008b906SSatish Balay 
1322682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat);
13235615d1e5SSatish Balay #undef __FUNC__
13245eea60f9SBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ" /* ADIC Ignore */
13258f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A)
1326682d7d0cSBarry Smith {
1327682d7d0cSBarry Smith   Mat_MPIAIJ *a   = (Mat_MPIAIJ*) A->data;
1328682d7d0cSBarry Smith 
1329682d7d0cSBarry Smith   if (!a->rank) return MatPrintHelp_SeqAIJ(a->A);
1330682d7d0cSBarry Smith   else return 0;
1331682d7d0cSBarry Smith }
1332682d7d0cSBarry Smith 
13335615d1e5SSatish Balay #undef __FUNC__
13345eea60f9SBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ" /* ADIC Ignore */
13358f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs)
13365a838052SSatish Balay {
13375a838052SSatish Balay   *bs = 1;
13385a838052SSatish Balay   return 0;
13395a838052SSatish Balay }
13405615d1e5SSatish Balay #undef __FUNC__
13415eea60f9SBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ" /* ADIC Ignore */
13428f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A)
1343bb5a7306SBarry Smith {
1344bb5a7306SBarry Smith   Mat_MPIAIJ *a   = (Mat_MPIAIJ*) A->data;
1345bb5a7306SBarry Smith   int        ierr;
1346bb5a7306SBarry Smith   ierr = MatSetUnfactored(a->A); CHKERRQ(ierr);
1347bb5a7306SBarry Smith   return 0;
1348bb5a7306SBarry Smith }
1349bb5a7306SBarry Smith 
13508f6be9afSLois Curfman McInnes extern int MatConvertSameType_MPIAIJ(Mat,Mat *,int);
13512f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int);
13520a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring);
13530a198c4cSBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **);
13548a729477SBarry Smith /* -------------------------------------------------------------------*/
13552ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ,
135639e00950SLois Curfman McInnes        MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ,
135744a69424SLois Curfman McInnes        MatMult_MPIAIJ,MatMultAdd_MPIAIJ,
135844a69424SLois Curfman McInnes        MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ,
135936ce4990SBarry Smith        0,0,
136036ce4990SBarry Smith        0,0,
136136ce4990SBarry Smith        0,0,
136244a69424SLois Curfman McInnes        MatRelax_MPIAIJ,
1363b7c46309SBarry Smith        MatTranspose_MPIAIJ,
1364a66be287SLois Curfman McInnes        MatGetInfo_MPIAIJ,0,
1365a008b906SSatish Balay        MatGetDiagonal_MPIAIJ,MatDiagonalScale_MPIAIJ,MatNorm_MPIAIJ,
1366ee50ffe9SBarry Smith        MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ,
13671eb62cbbSBarry Smith        0,
1368299609e3SLois Curfman McInnes        MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ,
136936ce4990SBarry Smith        0,0,0,0,
1370d1710a03SLois Curfman McInnes        MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ,
137136ce4990SBarry Smith        0,0,
137294a9d846SBarry Smith        0,0,MatConvertSameType_MPIAIJ,0,0,
1373b49de8d1SLois Curfman McInnes        0,0,0,
1374598137ffSSatish Balay        MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0,
1375052efed2SBarry Smith        MatPrintHelp_MPIAIJ,
13763b2fbd54SBarry Smith        MatScale_MPIAIJ,0,0,0,
13770a198c4cSBarry Smith        MatGetBlockSize_MPIAIJ,0,0,0,0,
1378bb5a7306SBarry Smith        MatFDColoringCreate_MPIAIJ,0,MatSetUnfactored_MPIAIJ};
137936ce4990SBarry Smith 
13808a729477SBarry Smith 
13815615d1e5SSatish Balay #undef __FUNC__
13825615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ"
13831987afe7SBarry Smith /*@C
1384ff756334SLois Curfman McInnes    MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format
13853a511b96SLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
13863a511b96SLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameters
13873a511b96SLois Curfman McInnes    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
13883a511b96SLois Curfman McInnes    performance can be increased by more than a factor of 50.
13898a729477SBarry Smith 
13908a729477SBarry Smith    Input Parameters:
13911eb62cbbSBarry Smith .  comm - MPI communicator
13927d3e4905SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
139392e8d321SLois Curfman McInnes            This value should be the same as the local size used in creating the
139492e8d321SLois Curfman McInnes            y vector for the matrix-vector product y = Ax.
139592e8d321SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
139692e8d321SLois Curfman McInnes            This value should be the same as the local size used in creating the
139792e8d321SLois Curfman McInnes            x vector for the matrix-vector product y = Ax.
13987d3e4905SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
139992e8d321SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1400ab693e5aSLois Curfman McInnes .  d_nz - number of nonzeros per row in diagonal portion of local submatrix
1401ff756334SLois Curfman McInnes            (same for all local rows)
14022bd5e0b2SLois Curfman McInnes .  d_nzz - array containing the number of nonzeros in the various rows of the
140392e8d321SLois Curfman McInnes            diagonal portion of the local submatrix (possibly different for each row)
14042bd5e0b2SLois Curfman McInnes            or PETSC_NULL. You must leave room for the diagonal entry even if
14052bd5e0b2SLois Curfman McInnes            it is zero.
14062bd5e0b2SLois Curfman McInnes .  o_nz - number of nonzeros per row in the off-diagonal portion of local
1407ab693e5aSLois Curfman McInnes            submatrix (same for all local rows).
14082bd5e0b2SLois Curfman McInnes .  o_nzz - array containing the number of nonzeros in the various rows of the
14092bd5e0b2SLois Curfman McInnes            off-diagonal portion of the local submatrix (possibly different for
14102bd5e0b2SLois Curfman McInnes            each row) or PETSC_NULL.
14118a729477SBarry Smith 
1412ff756334SLois Curfman McInnes    Output Parameter:
141344cd7ae7SLois Curfman McInnes .  A - the matrix
14148a729477SBarry Smith 
1415ff756334SLois Curfman McInnes    Notes:
1416ff756334SLois Curfman McInnes    The AIJ format (also called the Yale sparse matrix format or
1417ff756334SLois Curfman McInnes    compressed row storage), is fully compatible with standard Fortran 77
14180002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
14190002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
1420ff756334SLois Curfman McInnes 
1421ff756334SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
1422ff756334SLois Curfman McInnes    (possibly both).
1423ff756334SLois Curfman McInnes 
14245511cfe3SLois Curfman McInnes    By default, this format uses inodes (identical nodes) when possible.
14255511cfe3SLois Curfman McInnes    We search for consecutive rows with the same nonzero structure, thereby
14265511cfe3SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
14275511cfe3SLois Curfman McInnes 
14285511cfe3SLois Curfman McInnes    Options Database Keys:
14295511cfe3SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
14306e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
14316e62573dSLois Curfman McInnes $        (max limit=5)
14326e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
14336e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
14346e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
14355511cfe3SLois Curfman McInnes 
1436e0245417SLois Curfman McInnes    Storage Information:
1437e0245417SLois Curfman McInnes    For a square global matrix we define each processor's diagonal portion
1438e0245417SLois Curfman McInnes    to be its local rows and the corresponding columns (a square submatrix);
1439e0245417SLois Curfman McInnes    each processor's off-diagonal portion encompasses the remainder of the
1440e0245417SLois Curfman McInnes    local matrix (a rectangular submatrix).
1441e0245417SLois Curfman McInnes 
1442e0245417SLois Curfman McInnes    The user can specify preallocated storage for the diagonal part of
14435ace5be8SLois Curfman McInnes    the local submatrix with either d_nz or d_nnz (not both).  Set
14445ace5be8SLois Curfman McInnes    d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
14455ace5be8SLois Curfman McInnes    memory allocation.  Likewise, specify preallocated storage for the
14465ace5be8SLois Curfman McInnes    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
1447ff756334SLois Curfman McInnes 
14485511cfe3SLois Curfman McInnes    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
14495511cfe3SLois Curfman McInnes    the figure below we depict these three local rows and all columns (0-11).
14502191d07cSBarry Smith 
1451b810aeb4SBarry Smith $          0 1 2 3 4 5 6 7 8 9 10 11
1452b810aeb4SBarry Smith $         -------------------
14535511cfe3SLois Curfman McInnes $  row 3  |  o o o d d d o o o o o o
14545511cfe3SLois Curfman McInnes $  row 4  |  o o o d d d o o o o o o
14555511cfe3SLois Curfman McInnes $  row 5  |  o o o d d d o o o o o o
14565511cfe3SLois Curfman McInnes $         -------------------
1457b810aeb4SBarry Smith $
1458b810aeb4SBarry Smith 
14595511cfe3SLois Curfman McInnes    Thus, any entries in the d locations are stored in the d (diagonal)
14605511cfe3SLois Curfman McInnes    submatrix, and any entries in the o locations are stored in the
14615511cfe3SLois Curfman McInnes    o (off-diagonal) submatrix.  Note that the d and the o submatrices are
14625511cfe3SLois Curfman McInnes    stored simply in the MATSEQAIJ format for compressed row storage.
14635511cfe3SLois Curfman McInnes 
14645511cfe3SLois Curfman McInnes    Now d_nz should indicate the number of nonzeros per row in the d matrix,
14655511cfe3SLois Curfman McInnes    and o_nz should indicate the number of nonzeros per row in the o matrix.
14665511cfe3SLois Curfman McInnes    In general, for PDE problems in which most nonzeros are near the diagonal,
14673d323bbdSBarry Smith    one expects d_nz >> o_nz. For large problems you MUST preallocate memory
146892e8d321SLois Curfman McInnes    or you will get TERRIBLE performance; see the users' manual chapter on
14696da5968aSLois Curfman McInnes    matrices.
14703a511b96SLois Curfman McInnes 
1471dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel
1472ff756334SLois Curfman McInnes 
1473fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
14748a729477SBarry Smith @*/
14751eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,
147644cd7ae7SLois Curfman McInnes                     int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A)
14778a729477SBarry Smith {
147844cd7ae7SLois Curfman McInnes   Mat          B;
147944cd7ae7SLois Curfman McInnes   Mat_MPIAIJ   *b;
148036ce4990SBarry Smith   int          ierr, i,sum[2],work[2],size;
1481416022c9SBarry Smith 
148244cd7ae7SLois Curfman McInnes   *A = 0;
148344cd7ae7SLois Curfman McInnes   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATMPIAIJ,comm);
148444cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
148544cd7ae7SLois Curfman McInnes   B->data       = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b);
148644cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_MPIAIJ));
148744cd7ae7SLois Curfman McInnes   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
148836ce4990SBarry Smith   MPI_Comm_size(comm,&size);
148936ce4990SBarry Smith   if (size == 1) {
149036ce4990SBarry Smith     B->ops.getrowij          = MatGetRowIJ_MPIAIJ;
149136ce4990SBarry Smith     B->ops.restorerowij      = MatRestoreRowIJ_MPIAIJ;
149236ce4990SBarry Smith     B->ops.lufactorsymbolic  = MatLUFactorSymbolic_MPIAIJ;
149336ce4990SBarry Smith     B->ops.lufactornumeric   = MatLUFactorNumeric_MPIAIJ;
149436ce4990SBarry Smith     B->ops.lufactor          = MatLUFactor_MPIAIJ;
149536ce4990SBarry Smith     B->ops.solve             = MatSolve_MPIAIJ;
149636ce4990SBarry Smith     B->ops.solveadd          = MatSolveAdd_MPIAIJ;
149736ce4990SBarry Smith     B->ops.solvetrans        = MatSolveTrans_MPIAIJ;
149836ce4990SBarry Smith     B->ops.solvetransadd     = MatSolveTransAdd_MPIAIJ;
149936ce4990SBarry Smith     B->ops.ilufactorsymbolic = MatILUFactorSymbolic_MPIAIJ;
150036ce4990SBarry Smith   }
150144cd7ae7SLois Curfman McInnes   B->destroy    = MatDestroy_MPIAIJ;
150244cd7ae7SLois Curfman McInnes   B->view       = MatView_MPIAIJ;
150344cd7ae7SLois Curfman McInnes   B->factor     = 0;
150444cd7ae7SLois Curfman McInnes   B->assembled  = PETSC_FALSE;
150590f02eecSBarry Smith   B->mapping    = 0;
1506d6dfbf8fSBarry Smith 
150747794344SBarry Smith   B->insertmode = NOT_SET_VALUES;
150844cd7ae7SLois Curfman McInnes   MPI_Comm_rank(comm,&b->rank);
150944cd7ae7SLois Curfman McInnes   MPI_Comm_size(comm,&b->size);
15101eb62cbbSBarry Smith 
1511b4fd4287SBarry Smith   if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL))
1512e3372554SBarry Smith     SETERRQ(1,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz");
15131987afe7SBarry Smith 
1514dbd7a890SLois Curfman McInnes   if (M == PETSC_DECIDE || N == PETSC_DECIDE) {
15151eb62cbbSBarry Smith     work[0] = m; work[1] = n;
1516d65a2f8fSBarry Smith     MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm );
1517dbd7a890SLois Curfman McInnes     if (M == PETSC_DECIDE) M = sum[0];
1518dbd7a890SLois Curfman McInnes     if (N == PETSC_DECIDE) N = sum[1];
15191eb62cbbSBarry Smith   }
152044cd7ae7SLois Curfman McInnes   if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);}
152144cd7ae7SLois Curfman McInnes   if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);}
152244cd7ae7SLois Curfman McInnes   b->m = m; B->m = m;
152344cd7ae7SLois Curfman McInnes   b->n = n; B->n = n;
152444cd7ae7SLois Curfman McInnes   b->N = N; B->N = N;
152544cd7ae7SLois Curfman McInnes   b->M = M; B->M = M;
15261eb62cbbSBarry Smith 
15271eb62cbbSBarry Smith   /* build local table of row and column ownerships */
152844cd7ae7SLois Curfman McInnes   b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners);
152944cd7ae7SLois Curfman McInnes   PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ));
1530603f58a4SSatish Balay   b->cowners = b->rowners + b->size + 2;
153144cd7ae7SLois Curfman McInnes   MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);
153244cd7ae7SLois Curfman McInnes   b->rowners[0] = 0;
153344cd7ae7SLois Curfman McInnes   for ( i=2; i<=b->size; i++ ) {
153444cd7ae7SLois Curfman McInnes     b->rowners[i] += b->rowners[i-1];
15358a729477SBarry Smith   }
153644cd7ae7SLois Curfman McInnes   b->rstart = b->rowners[b->rank];
153744cd7ae7SLois Curfman McInnes   b->rend   = b->rowners[b->rank+1];
153844cd7ae7SLois Curfman McInnes   MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm);
153944cd7ae7SLois Curfman McInnes   b->cowners[0] = 0;
154044cd7ae7SLois Curfman McInnes   for ( i=2; i<=b->size; i++ ) {
154144cd7ae7SLois Curfman McInnes     b->cowners[i] += b->cowners[i-1];
15428a729477SBarry Smith   }
154344cd7ae7SLois Curfman McInnes   b->cstart = b->cowners[b->rank];
154444cd7ae7SLois Curfman McInnes   b->cend   = b->cowners[b->rank+1];
15458a729477SBarry Smith 
15465ace5be8SLois Curfman McInnes   if (d_nz == PETSC_DEFAULT) d_nz = 5;
154744cd7ae7SLois Curfman McInnes   ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr);
154844cd7ae7SLois Curfman McInnes   PLogObjectParent(B,b->A);
15497b8455f0SLois Curfman McInnes   if (o_nz == PETSC_DEFAULT) o_nz = 0;
155044cd7ae7SLois Curfman McInnes   ierr = MatCreateSeqAIJ(MPI_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr);
155144cd7ae7SLois Curfman McInnes   PLogObjectParent(B,b->B);
15528a729477SBarry Smith 
15531eb62cbbSBarry Smith   /* build cache for off array entries formed */
155444cd7ae7SLois Curfman McInnes   ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr);
155590f02eecSBarry Smith   b->donotstash  = 0;
155644cd7ae7SLois Curfman McInnes   b->colmap      = 0;
155744cd7ae7SLois Curfman McInnes   b->garray      = 0;
155844cd7ae7SLois Curfman McInnes   b->roworiented = 1;
15598a729477SBarry Smith 
15601eb62cbbSBarry Smith   /* stuff used for matrix vector multiply */
156144cd7ae7SLois Curfman McInnes   b->lvec      = 0;
156244cd7ae7SLois Curfman McInnes   b->Mvctx     = 0;
15638a729477SBarry Smith 
15647a0afa10SBarry Smith   /* stuff for MatGetRow() */
156544cd7ae7SLois Curfman McInnes   b->rowindices   = 0;
156644cd7ae7SLois Curfman McInnes   b->rowvalues    = 0;
156744cd7ae7SLois Curfman McInnes   b->getrowactive = PETSC_FALSE;
15687a0afa10SBarry Smith 
156944cd7ae7SLois Curfman McInnes   *A = B;
1570d6dfbf8fSBarry Smith   return 0;
1571d6dfbf8fSBarry Smith }
1572c74985f6SBarry Smith 
15735615d1e5SSatish Balay #undef __FUNC__
15745615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIAIJ"
15758f6be9afSLois Curfman McInnes int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues)
1576d6dfbf8fSBarry Smith {
1577d6dfbf8fSBarry Smith   Mat        mat;
1578416022c9SBarry Smith   Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data;
1579a1b97e82SLois Curfman McInnes   int        ierr, len=0, flg;
1580d6dfbf8fSBarry Smith 
1581416022c9SBarry Smith   *newmat       = 0;
15820452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm);
1583a5a9c739SBarry Smith   PLogObjectCreate(mat);
15840452661fSBarry Smith   mat->data       = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a);
1585416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
158644a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
158744a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
1588d6dfbf8fSBarry Smith   mat->factor     = matin->factor;
1589c456f294SBarry Smith   mat->assembled  = PETSC_TRUE;
1590d6dfbf8fSBarry Smith 
159144cd7ae7SLois Curfman McInnes   a->m = mat->m   = oldmat->m;
159244cd7ae7SLois Curfman McInnes   a->n = mat->n   = oldmat->n;
159344cd7ae7SLois Curfman McInnes   a->M = mat->M   = oldmat->M;
159444cd7ae7SLois Curfman McInnes   a->N = mat->N   = oldmat->N;
1595d6dfbf8fSBarry Smith 
1596416022c9SBarry Smith   a->rstart       = oldmat->rstart;
1597416022c9SBarry Smith   a->rend         = oldmat->rend;
1598416022c9SBarry Smith   a->cstart       = oldmat->cstart;
1599416022c9SBarry Smith   a->cend         = oldmat->cend;
160017699dbbSLois Curfman McInnes   a->size         = oldmat->size;
160117699dbbSLois Curfman McInnes   a->rank         = oldmat->rank;
160247794344SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1603bcd2baecSBarry Smith   a->rowvalues    = 0;
1604bcd2baecSBarry Smith   a->getrowactive = PETSC_FALSE;
1605d6dfbf8fSBarry Smith 
1606603f58a4SSatish Balay   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners);
1607603f58a4SSatish Balay   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ));
1608603f58a4SSatish Balay   a->cowners = a->rowners + a->size + 2;
1609603f58a4SSatish Balay   PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int));
1610416022c9SBarry Smith   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
16112ee70a88SLois Curfman McInnes   if (oldmat->colmap) {
16120452661fSBarry Smith     a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap);
1613416022c9SBarry Smith     PLogObjectMemory(mat,(a->N)*sizeof(int));
1614416022c9SBarry Smith     PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int));
1615416022c9SBarry Smith   } else a->colmap = 0;
16163f41c07dSBarry Smith   if (oldmat->garray) {
16173f41c07dSBarry Smith     len = ((Mat_SeqAIJ *) (oldmat->B->data))->n;
16183f41c07dSBarry Smith     a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray);
1619464493b3SBarry Smith     PLogObjectMemory(mat,len*sizeof(int));
16203f41c07dSBarry Smith     if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int));
1621416022c9SBarry Smith   } else a->garray = 0;
1622d6dfbf8fSBarry Smith 
1623416022c9SBarry Smith   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
1624416022c9SBarry Smith   PLogObjectParent(mat,a->lvec);
1625a56f8943SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
1626416022c9SBarry Smith   PLogObjectParent(mat,a->Mvctx);
1627416022c9SBarry Smith   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
1628416022c9SBarry Smith   PLogObjectParent(mat,a->A);
1629416022c9SBarry Smith   ierr =  MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr);
1630416022c9SBarry Smith   PLogObjectParent(mat,a->B);
16315dd7a6c7SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
163225cdf11fSBarry Smith   if (flg) {
1633682d7d0cSBarry Smith     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
1634682d7d0cSBarry Smith   }
16358a729477SBarry Smith   *newmat = mat;
16368a729477SBarry Smith   return 0;
16378a729477SBarry Smith }
1638416022c9SBarry Smith 
163977c4ece6SBarry Smith #include "sys.h"
1640416022c9SBarry Smith 
16415615d1e5SSatish Balay #undef __FUNC__
16425615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ"
164319bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat)
1644416022c9SBarry Smith {
1645d65a2f8fSBarry Smith   Mat          A;
1646416022c9SBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
1647d65a2f8fSBarry Smith   Scalar       *vals,*svals;
164819bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1649416022c9SBarry Smith   MPI_Status   status;
165017699dbbSLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
1651d65a2f8fSBarry Smith   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
165219bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
1653416022c9SBarry Smith 
165417699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank);
165517699dbbSLois Curfman McInnes   if (!rank) {
165690ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
165777c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr);
1658e3372554SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object");
1659416022c9SBarry Smith   }
1660416022c9SBarry Smith 
1661416022c9SBarry Smith   MPI_Bcast(header+1,3,MPI_INT,0,comm);
1662416022c9SBarry Smith   M = header[1]; N = header[2];
1663416022c9SBarry Smith   /* determine ownership of all rows */
166417699dbbSLois Curfman McInnes   m = M/size + ((M % size) > rank);
16650452661fSBarry Smith   rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
1666416022c9SBarry Smith   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1667416022c9SBarry Smith   rowners[0] = 0;
166817699dbbSLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
1669416022c9SBarry Smith     rowners[i] += rowners[i-1];
1670416022c9SBarry Smith   }
167117699dbbSLois Curfman McInnes   rstart = rowners[rank];
167217699dbbSLois Curfman McInnes   rend   = rowners[rank+1];
1673416022c9SBarry Smith 
1674416022c9SBarry Smith   /* distribute row lengths to all processors */
16750452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
1676416022c9SBarry Smith   offlens = ourlens + (rend-rstart);
167717699dbbSLois Curfman McInnes   if (!rank) {
16780452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
167977c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
16800452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts);
168117699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1682d65a2f8fSBarry Smith     MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
16830452661fSBarry Smith     PetscFree(sndcounts);
1684416022c9SBarry Smith   }
1685416022c9SBarry Smith   else {
1686416022c9SBarry Smith     MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);
1687416022c9SBarry Smith   }
1688416022c9SBarry Smith 
168917699dbbSLois Curfman McInnes   if (!rank) {
1690416022c9SBarry Smith     /* calculate the number of nonzeros on each processor */
16910452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz);
1692cddf8d76SBarry Smith     PetscMemzero(procsnz,size*sizeof(int));
169317699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) {
1694416022c9SBarry Smith       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
1695416022c9SBarry Smith         procsnz[i] += rowlengths[j];
1696416022c9SBarry Smith       }
1697416022c9SBarry Smith     }
16980452661fSBarry Smith     PetscFree(rowlengths);
1699416022c9SBarry Smith 
1700416022c9SBarry Smith     /* determine max buffer needed and allocate it */
1701416022c9SBarry Smith     maxnz = 0;
170217699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) {
17030452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
1704416022c9SBarry Smith     }
17050452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols);
1706416022c9SBarry Smith 
1707416022c9SBarry Smith     /* read in my part of the matrix column indices  */
1708416022c9SBarry Smith     nz = procsnz[0];
17090452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
171077c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr);
1711d65a2f8fSBarry Smith 
1712d65a2f8fSBarry Smith     /* read in every one elses and ship off */
171317699dbbSLois Curfman McInnes     for ( i=1; i<size; i++ ) {
1714d65a2f8fSBarry Smith       nz = procsnz[i];
171577c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
1716d65a2f8fSBarry Smith       MPI_Send(cols,nz,MPI_INT,i,tag,comm);
1717d65a2f8fSBarry Smith     }
17180452661fSBarry Smith     PetscFree(cols);
1719416022c9SBarry Smith   }
1720416022c9SBarry Smith   else {
1721416022c9SBarry Smith     /* determine buffer space needed for message */
1722416022c9SBarry Smith     nz = 0;
1723416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
1724416022c9SBarry Smith       nz += ourlens[i];
1725416022c9SBarry Smith     }
17260452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
1727416022c9SBarry Smith 
1728416022c9SBarry Smith     /* receive message of column indices*/
1729d65a2f8fSBarry Smith     MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);
1730416022c9SBarry Smith     MPI_Get_count(&status,MPI_INT,&maxnz);
1731e3372554SBarry Smith     if (maxnz != nz) SETERRQ(1,0,"something is wrong with file");
1732416022c9SBarry Smith   }
1733416022c9SBarry Smith 
1734416022c9SBarry Smith   /* loop over local rows, determining number of off diagonal entries */
1735cddf8d76SBarry Smith   PetscMemzero(offlens,m*sizeof(int));
1736416022c9SBarry Smith   jj = 0;
1737416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1738416022c9SBarry Smith     for ( j=0; j<ourlens[i]; j++ ) {
1739d65a2f8fSBarry Smith       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
1740416022c9SBarry Smith       jj++;
1741416022c9SBarry Smith     }
1742416022c9SBarry Smith   }
1743d65a2f8fSBarry Smith 
1744d65a2f8fSBarry Smith   /* create our matrix */
1745416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1746416022c9SBarry Smith     ourlens[i] -= offlens[i];
1747416022c9SBarry Smith   }
1748d65a2f8fSBarry Smith   ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr);
1749d65a2f8fSBarry Smith   A = *newmat;
17506d4a8577SBarry Smith   MatSetOption(A,MAT_COLUMNS_SORTED);
1751d65a2f8fSBarry Smith   for ( i=0; i<m; i++ ) {
1752d65a2f8fSBarry Smith     ourlens[i] += offlens[i];
1753d65a2f8fSBarry Smith   }
1754416022c9SBarry Smith 
175517699dbbSLois Curfman McInnes   if (!rank) {
17560452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
1757416022c9SBarry Smith 
1758416022c9SBarry Smith     /* read in my part of the matrix numerical values  */
1759416022c9SBarry Smith     nz = procsnz[0];
176077c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
1761d65a2f8fSBarry Smith 
1762d65a2f8fSBarry Smith     /* insert into matrix */
1763d65a2f8fSBarry Smith     jj      = rstart;
1764d65a2f8fSBarry Smith     smycols = mycols;
1765d65a2f8fSBarry Smith     svals   = vals;
1766d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1767d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1768d65a2f8fSBarry Smith       smycols += ourlens[i];
1769d65a2f8fSBarry Smith       svals   += ourlens[i];
1770d65a2f8fSBarry Smith       jj++;
1771416022c9SBarry Smith     }
1772416022c9SBarry Smith 
1773d65a2f8fSBarry Smith     /* read in other processors and ship out */
177417699dbbSLois Curfman McInnes     for ( i=1; i<size; i++ ) {
1775416022c9SBarry Smith       nz = procsnz[i];
177677c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
1777d65a2f8fSBarry Smith       MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
1778416022c9SBarry Smith     }
17790452661fSBarry Smith     PetscFree(procsnz);
1780416022c9SBarry Smith   }
1781d65a2f8fSBarry Smith   else {
1782d65a2f8fSBarry Smith     /* receive numeric values */
17830452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
1784416022c9SBarry Smith 
1785d65a2f8fSBarry Smith     /* receive message of values*/
1786d65a2f8fSBarry Smith     MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
1787d65a2f8fSBarry Smith     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1788e3372554SBarry Smith     if (maxnz != nz) SETERRQ(1,0,"something is wrong with file");
1789d65a2f8fSBarry Smith 
1790d65a2f8fSBarry Smith     /* insert into matrix */
1791d65a2f8fSBarry Smith     jj      = rstart;
1792d65a2f8fSBarry Smith     smycols = mycols;
1793d65a2f8fSBarry Smith     svals   = vals;
1794d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1795d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1796d65a2f8fSBarry Smith       smycols += ourlens[i];
1797d65a2f8fSBarry Smith       svals   += ourlens[i];
1798d65a2f8fSBarry Smith       jj++;
1799d65a2f8fSBarry Smith     }
1800d65a2f8fSBarry Smith   }
18010452661fSBarry Smith   PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners);
1802d65a2f8fSBarry Smith 
18036d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
18046d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1805416022c9SBarry Smith   return 0;
1806416022c9SBarry Smith }
1807