xref: /petsc/src/mat/impls/aij/mpi/mpiaij.c (revision b5bd3ad58742d10278d184b667b60ae1619f8109)
1b57c8cebSBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*b5bd3ad5SBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.244 1998/05/11 21:53:56 bsmith Exp bsmith $";
4cb512458SBarry Smith #endif
58a729477SBarry Smith 
63369ce9aSBarry Smith #include "pinclude/pviewer.h"
770f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h"
8f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
9d9942c19SSatish Balay #include "src/inline/spops.h"
108a729477SBarry Smith 
119e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column
129e25ed09SBarry Smith number to the local number in the off-diagonal part of the local
139e25ed09SBarry Smith storage of the matrix.  This is done in a non scable way since the
149e25ed09SBarry Smith length of colmap equals the global matrix length.
159e25ed09SBarry Smith */
165615d1e5SSatish Balay #undef __FUNC__
17d4bb536fSBarry Smith #define __FUNC__ "CreateColmap_MPIAIJ_Private"
180a198c4cSBarry Smith int CreateColmap_MPIAIJ_Private(Mat mat)
199e25ed09SBarry Smith {
2044a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
21ec8511deSBarry Smith   Mat_SeqAIJ *B = (Mat_SeqAIJ*) aij->B->data;
22905e6a2fSBarry Smith   int        n = B->n,i;
23dbb450caSBarry Smith 
243a40ed3dSBarry Smith   PetscFunctionBegin;
25758f045eSSatish Balay   aij->colmap = (int *) PetscMalloc((aij->N+1)*sizeof(int));CHKPTRQ(aij->colmap);
26464493b3SBarry Smith   PLogObjectMemory(mat,aij->N*sizeof(int));
27cddf8d76SBarry Smith   PetscMemzero(aij->colmap,aij->N*sizeof(int));
28905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1;
293a40ed3dSBarry Smith   PetscFunctionReturn(0);
309e25ed09SBarry Smith }
319e25ed09SBarry Smith 
322493cbb0SBarry Smith extern int DisAssemble_MPIAIJ(Mat);
332493cbb0SBarry Smith 
340520107fSSatish Balay #define CHUNKSIZE   15
3530770e4dSSatish Balay #define MatSetValues_SeqAIJ_A_Private(row,col,value,addv) \
360520107fSSatish Balay { \
370520107fSSatish Balay  \
380520107fSSatish Balay     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift; \
3930770e4dSSatish Balay     rmax = aimax[row]; nrow = ailen[row];  \
40f5e9677aSSatish Balay     col1 = col - shift; \
41f5e9677aSSatish Balay      \
42ba4e3ef2SSatish Balay     low = 0; high = nrow; \
43ba4e3ef2SSatish Balay     while (high-low > 5) { \
44ba4e3ef2SSatish Balay       t = (low+high)/2; \
45ba4e3ef2SSatish Balay       if (rp[t] > col) high = t; \
46ba4e3ef2SSatish Balay       else             low  = t; \
47ba4e3ef2SSatish Balay     } \
480520107fSSatish Balay       for ( _i=0; _i<nrow; _i++ ) { \
49f5e9677aSSatish Balay         if (rp[_i] > col1) break; \
50f5e9677aSSatish Balay         if (rp[_i] == col1) { \
510520107fSSatish Balay           if (addv == ADD_VALUES) ap[_i] += value;   \
520520107fSSatish Balay           else                  ap[_i] = value; \
5330770e4dSSatish Balay           goto a_noinsert; \
540520107fSSatish Balay         } \
550520107fSSatish Balay       }  \
5689280ab3SLois Curfman McInnes       if (nonew == 1) goto a_noinsert; \
57a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero into matrix"); \
580520107fSSatish Balay       if (nrow >= rmax) { \
590520107fSSatish Balay         /* there is no extra room in row, therefore enlarge */ \
600520107fSSatish Balay         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; \
610520107fSSatish Balay         Scalar *new_a; \
620520107fSSatish Balay  \
63a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); \
6489280ab3SLois Curfman McInnes  \
650520107fSSatish Balay         /* malloc new storage space */ \
660520107fSSatish Balay         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); \
670520107fSSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \
680520107fSSatish Balay         new_j   = (int *) (new_a + new_nz); \
690520107fSSatish Balay         new_i   = new_j + new_nz; \
700520107fSSatish Balay  \
710520107fSSatish Balay         /* copy over old data into new slots */ \
720520107fSSatish Balay         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} \
730520107fSSatish Balay         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} \
740520107fSSatish Balay         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); \
750520107fSSatish Balay         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); \
760520107fSSatish Balay         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, \
770520107fSSatish Balay                                                            len*sizeof(int)); \
780520107fSSatish Balay         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); \
790520107fSSatish Balay         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, \
800520107fSSatish Balay                                                            len*sizeof(Scalar));  \
810520107fSSatish Balay         /* free up old matrix storage */ \
82f5e9677aSSatish Balay  \
830520107fSSatish Balay         PetscFree(a->a);  \
840520107fSSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} \
850520107fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;  \
860520107fSSatish Balay         a->singlemalloc = 1; \
870520107fSSatish Balay  \
880520107fSSatish Balay         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift; \
8930770e4dSSatish Balay         rmax = aimax[row] = aimax[row] + CHUNKSIZE; \
900520107fSSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \
910520107fSSatish Balay         a->maxnz += CHUNKSIZE; \
920520107fSSatish Balay         a->reallocs++; \
930520107fSSatish Balay       } \
940520107fSSatish Balay       N = nrow++ - 1; a->nz++; \
950520107fSSatish Balay       /* shift up all the later entries in this row */ \
960520107fSSatish Balay       for ( ii=N; ii>=_i; ii-- ) { \
970520107fSSatish Balay         rp[ii+1] = rp[ii]; \
980520107fSSatish Balay         ap[ii+1] = ap[ii]; \
990520107fSSatish Balay       } \
100f5e9677aSSatish Balay       rp[_i] = col1;  \
1010520107fSSatish Balay       ap[_i] = value;  \
10230770e4dSSatish Balay       a_noinsert: ; \
1030520107fSSatish Balay       ailen[row] = nrow; \
1040520107fSSatish Balay }
1050a198c4cSBarry Smith 
10630770e4dSSatish Balay #define MatSetValues_SeqAIJ_B_Private(row,col,value,addv) \
10730770e4dSSatish Balay { \
10830770e4dSSatish Balay  \
10930770e4dSSatish Balay     rp   = bj + bi[row] + shift; ap = ba + bi[row] + shift; \
11030770e4dSSatish Balay     rmax = bimax[row]; nrow = bilen[row];  \
11130770e4dSSatish Balay     col1 = col - shift; \
11230770e4dSSatish Balay      \
113ba4e3ef2SSatish Balay     low = 0; high = nrow; \
114ba4e3ef2SSatish Balay     while (high-low > 5) { \
115ba4e3ef2SSatish Balay       t = (low+high)/2; \
116ba4e3ef2SSatish Balay       if (rp[t] > col) high = t; \
117ba4e3ef2SSatish Balay       else             low  = t; \
118ba4e3ef2SSatish Balay     } \
11930770e4dSSatish Balay        for ( _i=0; _i<nrow; _i++ ) { \
12030770e4dSSatish Balay         if (rp[_i] > col1) break; \
12130770e4dSSatish Balay         if (rp[_i] == col1) { \
12230770e4dSSatish Balay           if (addv == ADD_VALUES) ap[_i] += value;   \
12330770e4dSSatish Balay           else                  ap[_i] = value; \
12430770e4dSSatish Balay           goto b_noinsert; \
12530770e4dSSatish Balay         } \
12630770e4dSSatish Balay       }  \
12789280ab3SLois Curfman McInnes       if (nonew == 1) goto b_noinsert; \
128a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero into matrix"); \
12930770e4dSSatish Balay       if (nrow >= rmax) { \
13030770e4dSSatish Balay         /* there is no extra room in row, therefore enlarge */ \
13174c639caSSatish Balay         int    new_nz = bi[b->m] + CHUNKSIZE,len,*new_i,*new_j; \
13230770e4dSSatish Balay         Scalar *new_a; \
13330770e4dSSatish Balay  \
134a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); \
13589280ab3SLois Curfman McInnes  \
13630770e4dSSatish Balay         /* malloc new storage space */ \
13774c639caSSatish Balay         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(b->m+1)*sizeof(int); \
13830770e4dSSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \
13930770e4dSSatish Balay         new_j   = (int *) (new_a + new_nz); \
14030770e4dSSatish Balay         new_i   = new_j + new_nz; \
14130770e4dSSatish Balay  \
14230770e4dSSatish Balay         /* copy over old data into new slots */ \
14330770e4dSSatish Balay         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = bi[ii];} \
14474c639caSSatish Balay         for ( ii=row+1; ii<b->m+1; ii++ ) {new_i[ii] = bi[ii]+CHUNKSIZE;} \
14530770e4dSSatish Balay         PetscMemcpy(new_j,bj,(bi[row]+nrow+shift)*sizeof(int)); \
14630770e4dSSatish Balay         len = (new_nz - CHUNKSIZE - bi[row] - nrow - shift); \
14730770e4dSSatish Balay         PetscMemcpy(new_j+bi[row]+shift+nrow+CHUNKSIZE,bj+bi[row]+shift+nrow, \
14830770e4dSSatish Balay                                                            len*sizeof(int)); \
14930770e4dSSatish Balay         PetscMemcpy(new_a,ba,(bi[row]+nrow+shift)*sizeof(Scalar)); \
15030770e4dSSatish Balay         PetscMemcpy(new_a+bi[row]+shift+nrow+CHUNKSIZE,ba+bi[row]+shift+nrow, \
15130770e4dSSatish Balay                                                            len*sizeof(Scalar));  \
15230770e4dSSatish Balay         /* free up old matrix storage */ \
15330770e4dSSatish Balay  \
15474c639caSSatish Balay         PetscFree(b->a);  \
15574c639caSSatish Balay         if (!b->singlemalloc) {PetscFree(b->i);PetscFree(b->j);} \
15674c639caSSatish Balay         ba = b->a = new_a; bi = b->i = new_i; bj = b->j = new_j;  \
15774c639caSSatish Balay         b->singlemalloc = 1; \
15830770e4dSSatish Balay  \
15930770e4dSSatish Balay         rp   = bj + bi[row] + shift; ap = ba + bi[row] + shift; \
16030770e4dSSatish Balay         rmax = bimax[row] = bimax[row] + CHUNKSIZE; \
16174c639caSSatish Balay         PLogObjectMemory(B,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \
16274c639caSSatish Balay         b->maxnz += CHUNKSIZE; \
16374c639caSSatish Balay         b->reallocs++; \
16430770e4dSSatish Balay       } \
16574c639caSSatish Balay       N = nrow++ - 1; b->nz++; \
16630770e4dSSatish Balay       /* shift up all the later entries in this row */ \
16730770e4dSSatish Balay       for ( ii=N; ii>=_i; ii-- ) { \
16830770e4dSSatish Balay         rp[ii+1] = rp[ii]; \
16930770e4dSSatish Balay         ap[ii+1] = ap[ii]; \
17030770e4dSSatish Balay       } \
17130770e4dSSatish Balay       rp[_i] = col1;  \
17230770e4dSSatish Balay       ap[_i] = value;  \
17330770e4dSSatish Balay       b_noinsert: ; \
17430770e4dSSatish Balay       bilen[row] = nrow; \
17530770e4dSSatish Balay }
17630770e4dSSatish Balay 
1770520107fSSatish Balay extern int MatSetValues_SeqAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode);
1785615d1e5SSatish Balay #undef __FUNC__
1795615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIAIJ"
1808f6be9afSLois Curfman McInnes int MatSetValues_MPIAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv)
1818a729477SBarry Smith {
18244a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1834b0e389bSBarry Smith   Scalar     value;
1841eb62cbbSBarry Smith   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
1851eb62cbbSBarry Smith   int        cstart = aij->cstart, cend = aij->cend,row,col;
186905e6a2fSBarry Smith   int        roworiented = aij->roworiented;
1878a729477SBarry Smith 
1880520107fSSatish Balay   /* Some Variables required in the macro */
1894ee7247eSSatish Balay   Mat        A = aij->A;
1904ee7247eSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
19130770e4dSSatish Balay   int        *aimax = a->imax, *ai = a->i, *ailen = a->ilen,*aj = a->j;
19230770e4dSSatish Balay   Scalar     *aa = a->a;
19330770e4dSSatish Balay 
19430770e4dSSatish Balay   Mat        B = aij->B;
19530770e4dSSatish Balay   Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data;
19630770e4dSSatish Balay   int        *bimax = b->imax, *bi = b->i, *bilen = b->ilen,*bj = b->j;
19730770e4dSSatish Balay   Scalar     *ba = b->a;
19830770e4dSSatish Balay 
199ba4e3ef2SSatish Balay   int        *rp,ii,nrow,_i,rmax, N, col1,low,high,t;
20030770e4dSSatish Balay   int        nonew = a->nonew,shift = a->indexshift;
20130770e4dSSatish Balay   Scalar     *ap;
2024ee7247eSSatish Balay 
2033a40ed3dSBarry Smith   PetscFunctionBegin;
2048a729477SBarry Smith   for ( i=0; i<m; i++ ) {
2053a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
206a8c6a408SBarry Smith     if (im[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
207a8c6a408SBarry Smith     if (im[i] >= aij->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
2080a198c4cSBarry Smith #endif
2094b0e389bSBarry Smith     if (im[i] >= rstart && im[i] < rend) {
2104b0e389bSBarry Smith       row = im[i] - rstart;
2111eb62cbbSBarry Smith       for ( j=0; j<n; j++ ) {
2124b0e389bSBarry Smith         if (in[j] >= cstart && in[j] < cend){
2134b0e389bSBarry Smith           col = in[j] - cstart;
2144b0e389bSBarry Smith           if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
21530770e4dSSatish Balay           MatSetValues_SeqAIJ_A_Private(row,col,value,addv);
2160520107fSSatish Balay           /* ierr = MatSetValues_SeqAIJ(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */
2171eb62cbbSBarry Smith         }
2183a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
219a8c6a408SBarry Smith         else if (in[j] < 0) {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");}
220a8c6a408SBarry Smith         else if (in[j] >= aij->N) {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");}
2210a198c4cSBarry Smith #endif
2221eb62cbbSBarry Smith         else {
223227d817aSBarry Smith           if (mat->was_assembled) {
224905e6a2fSBarry Smith             if (!aij->colmap) {
225905e6a2fSBarry Smith               ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
226905e6a2fSBarry Smith             }
227905e6a2fSBarry Smith             col = aij->colmap[in[j]] - 1;
228ec8511deSBarry Smith             if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) {
2292493cbb0SBarry Smith               ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
2304b0e389bSBarry Smith               col =  in[j];
2319bf004c3SSatish Balay               /* Reinitialize the variables required by MatSetValues_SeqAIJ_B_Private() */
232f9508a3cSSatish Balay               B = aij->B;
233f9508a3cSSatish Balay               b = (Mat_SeqAIJ *) B->data;
234f9508a3cSSatish Balay               bimax = b->imax; bi = b->i; bilen = b->ilen; bj = b->j;
235f9508a3cSSatish Balay               ba = b->a;
236d6dfbf8fSBarry Smith             }
237c48de900SBarry Smith           } else col = in[j];
2384b0e389bSBarry Smith           if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
23930770e4dSSatish Balay           MatSetValues_SeqAIJ_B_Private(row,col,value,addv);
24030770e4dSSatish Balay           /* ierr = MatSetValues_SeqAIJ(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */
2411eb62cbbSBarry Smith         }
2421eb62cbbSBarry Smith       }
2431eb62cbbSBarry Smith     }
2441eb62cbbSBarry Smith     else {
24590f02eecSBarry Smith       if (roworiented && !aij->donotstash) {
2464b0e389bSBarry Smith         ierr = StashValues_Private(&aij->stash,im[i],n,in,v+i*n,addv);CHKERRQ(ierr);
2474b0e389bSBarry Smith       }
2484b0e389bSBarry Smith       else {
24990f02eecSBarry Smith         if (!aij->donotstash) {
2504b0e389bSBarry Smith           row = im[i];
2514b0e389bSBarry Smith           for ( j=0; j<n; j++ ) {
2524b0e389bSBarry Smith             ierr = StashValues_Private(&aij->stash,row,1,in+j,v+i+j*m,addv);CHKERRQ(ierr);
2534b0e389bSBarry Smith           }
2544b0e389bSBarry Smith         }
2551eb62cbbSBarry Smith       }
2568a729477SBarry Smith     }
25790f02eecSBarry Smith   }
2583a40ed3dSBarry Smith   PetscFunctionReturn(0);
2598a729477SBarry Smith }
2608a729477SBarry Smith 
2615615d1e5SSatish Balay #undef __FUNC__
2625615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIAIJ"
2638f6be9afSLois Curfman McInnes int MatGetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
264b49de8d1SLois Curfman McInnes {
265b49de8d1SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
266b49de8d1SLois Curfman McInnes   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
267b49de8d1SLois Curfman McInnes   int        cstart = aij->cstart, cend = aij->cend,row,col;
268b49de8d1SLois Curfman McInnes 
2693a40ed3dSBarry Smith   PetscFunctionBegin;
270b49de8d1SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
271a8c6a408SBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
272a8c6a408SBarry Smith     if (idxm[i] >= aij->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
273b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
274b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
275b49de8d1SLois Curfman McInnes       for ( j=0; j<n; j++ ) {
276a8c6a408SBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
277a8c6a408SBarry Smith         if (idxn[j] >= aij->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
278b49de8d1SLois Curfman McInnes         if (idxn[j] >= cstart && idxn[j] < cend){
279b49de8d1SLois Curfman McInnes           col = idxn[j] - cstart;
280b49de8d1SLois Curfman McInnes           ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr);
281b49de8d1SLois Curfman McInnes         }
282b49de8d1SLois Curfman McInnes         else {
283905e6a2fSBarry Smith           if (!aij->colmap) {
284905e6a2fSBarry Smith             ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
285905e6a2fSBarry Smith           }
286905e6a2fSBarry Smith           col = aij->colmap[idxn[j]] - 1;
287e60e1c95SSatish Balay           if ((col < 0) || (aij->garray[col] != idxn[j])) *(v+i*n+j) = 0.0;
288d9d09a02SSatish Balay           else {
289b49de8d1SLois Curfman McInnes             ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr);
290b49de8d1SLois Curfman McInnes           }
291b49de8d1SLois Curfman McInnes         }
292b49de8d1SLois Curfman McInnes       }
293a8c6a408SBarry Smith     } else {
294a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported");
295b49de8d1SLois Curfman McInnes     }
296b49de8d1SLois Curfman McInnes   }
2973a40ed3dSBarry Smith   PetscFunctionReturn(0);
298b49de8d1SLois Curfman McInnes }
299b49de8d1SLois Curfman McInnes 
3005615d1e5SSatish Balay #undef __FUNC__
3015615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ"
3028f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode)
3038a729477SBarry Smith {
30444a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
305d6dfbf8fSBarry Smith   MPI_Comm    comm = mat->comm;
30617699dbbSLois Curfman McInnes   int         size = aij->size, *owners = aij->rowners;
30717699dbbSLois Curfman McInnes   int         rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr;
3081eb62cbbSBarry Smith   MPI_Request *send_waits,*recv_waits;
3096abc6512SBarry Smith   int         *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
3101eb62cbbSBarry Smith   InsertMode  addv;
3111eb62cbbSBarry Smith   Scalar      *rvalues,*svalues;
3121eb62cbbSBarry Smith 
3133a40ed3dSBarry Smith   PetscFunctionBegin;
314*b5bd3ad5SBarry Smith   if (aij->donotstash) {
315*b5bd3ad5SBarry Smith     aij->svalues    = 0; aij->rvalues    = 0;
316*b5bd3ad5SBarry Smith     aij->nsends     = 0; aij->nrecvs     = 0;
317*b5bd3ad5SBarry Smith     aij->send_waits = 0; aij->recv_waits = 0;
318*b5bd3ad5SBarry Smith     aij->rmax       = 0;
319*b5bd3ad5SBarry Smith     PetscFunctionReturn(0);
320*b5bd3ad5SBarry Smith   }
321*b5bd3ad5SBarry Smith 
3221eb62cbbSBarry Smith   /* make sure all processors are either in INSERTMODE or ADDMODE */
323ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
324dbb450caSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
325a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Some processors inserted others added");
3261eb62cbbSBarry Smith   }
32747794344SBarry Smith   mat->insertmode = addv; /* in case this processor had no cache */
3281eb62cbbSBarry Smith 
3291eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
3300452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
331cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
3320452661fSBarry Smith   owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
3331eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
3341eb62cbbSBarry Smith     idx = aij->stash.idx[i];
33517699dbbSLois Curfman McInnes     for ( j=0; j<size; j++ ) {
3361eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
3371eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
3388a729477SBarry Smith       }
3398a729477SBarry Smith     }
3408a729477SBarry Smith   }
34117699dbbSLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
3421eb62cbbSBarry Smith 
3431eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
3440452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
345ca161407SBarry Smith   ierr = MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
34617699dbbSLois Curfman McInnes   nreceives = work[rank];
347ca161407SBarry Smith   ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
34817699dbbSLois Curfman McInnes   nmax = work[rank];
3490452661fSBarry Smith   PetscFree(work);
3501eb62cbbSBarry Smith 
3511eb62cbbSBarry Smith   /* post receives:
3521eb62cbbSBarry Smith        1) each message will consist of ordered pairs
3531eb62cbbSBarry Smith      (global index,value) we store the global index as a double
354d6dfbf8fSBarry Smith      to simplify the message passing.
3551eb62cbbSBarry Smith        2) since we don't know how long each individual message is we
3561eb62cbbSBarry Smith      allocate the largest needed buffer for each receive. Potentially
3571eb62cbbSBarry Smith      this is a lot of wasted space.
3581eb62cbbSBarry Smith 
3591eb62cbbSBarry Smith 
3601eb62cbbSBarry Smith        This could be done better.
3611eb62cbbSBarry Smith   */
362ca161407SBarry Smith   rvalues    = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));CHKPTRQ(rvalues);
363ca161407SBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
3641eb62cbbSBarry Smith   for ( i=0; i<nreceives; i++ ) {
365ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,
366ca161407SBarry Smith               comm,recv_waits+i);CHKERRQ(ierr);
3671eb62cbbSBarry Smith   }
3681eb62cbbSBarry Smith 
3691eb62cbbSBarry Smith   /* do sends:
3701eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
3711eb62cbbSBarry Smith          the ith processor
3721eb62cbbSBarry Smith   */
3730452661fSBarry Smith   svalues    = (Scalar *) PetscMalloc(3*(aij->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues);
374ca161407SBarry Smith   send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
3750452661fSBarry Smith   starts     = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts);
3761eb62cbbSBarry Smith   starts[0]  = 0;
37717699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3781eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
3791eb62cbbSBarry Smith     svalues[3*starts[owner[i]]]       = (Scalar)  aij->stash.idx[i];
3801eb62cbbSBarry Smith     svalues[3*starts[owner[i]]+1]     = (Scalar)  aij->stash.idy[i];
3811eb62cbbSBarry Smith     svalues[3*(starts[owner[i]]++)+2] =  aij->stash.array[i];
3821eb62cbbSBarry Smith   }
3830452661fSBarry Smith   PetscFree(owner);
3841eb62cbbSBarry Smith   starts[0] = 0;
38517699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3861eb62cbbSBarry Smith   count = 0;
38717699dbbSLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3881eb62cbbSBarry Smith     if (procs[i]) {
389ca161407SBarry Smith       ierr = MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3901eb62cbbSBarry Smith     }
3911eb62cbbSBarry Smith   }
392*b5bd3ad5SBarry Smith   PetscFree(starts);
393*b5bd3ad5SBarry Smith   PetscFree(nprocs);
3941eb62cbbSBarry Smith 
3951eb62cbbSBarry Smith   /* Free cache space */
39610a665d1SBarry Smith   PLogInfo(aij->A,"MatAssemblyBegin_MPIAIJ:Number of off-processor values %d\n",aij->stash.n);
39778b31e54SBarry Smith   ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr);
3981eb62cbbSBarry Smith 
3991eb62cbbSBarry Smith   aij->svalues    = svalues;    aij->rvalues    = rvalues;
4001eb62cbbSBarry Smith   aij->nsends     = nsends;     aij->nrecvs     = nreceives;
4011eb62cbbSBarry Smith   aij->send_waits = send_waits; aij->recv_waits = recv_waits;
4021eb62cbbSBarry Smith   aij->rmax       = nmax;
4031eb62cbbSBarry Smith 
4043a40ed3dSBarry Smith   PetscFunctionReturn(0);
4051eb62cbbSBarry Smith }
40644a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat);
4071eb62cbbSBarry Smith 
4085615d1e5SSatish Balay #undef __FUNC__
4095615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIAIJ"
4108f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode)
4111eb62cbbSBarry Smith {
41244a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
4131eb62cbbSBarry Smith   MPI_Status  *send_status,recv_status;
414416022c9SBarry Smith   int         imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr;
415905e6a2fSBarry Smith   int         row,col,other_disassembled;
4161eb62cbbSBarry Smith   Scalar      *values,val;
41747794344SBarry Smith   InsertMode  addv = mat->insertmode;
4181eb62cbbSBarry Smith 
4193a40ed3dSBarry Smith   PetscFunctionBegin;
420*b5bd3ad5SBarry Smith 
4211eb62cbbSBarry Smith   /*  wait on receives */
4221eb62cbbSBarry Smith   while (count) {
423ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4241eb62cbbSBarry Smith     /* unpack receives into our local space */
425d6dfbf8fSBarry Smith     values = aij->rvalues + 3*imdex*aij->rmax;
426ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_SCALAR,&n);CHKERRQ(ierr);
4271eb62cbbSBarry Smith     n = n/3;
4281eb62cbbSBarry Smith     for ( i=0; i<n; i++ ) {
429227d817aSBarry Smith       row = (int) PetscReal(values[3*i]) - aij->rstart;
430227d817aSBarry Smith       col = (int) PetscReal(values[3*i+1]);
4311eb62cbbSBarry Smith       val = values[3*i+2];
4321eb62cbbSBarry Smith       if (col >= aij->cstart && col < aij->cend) {
4331eb62cbbSBarry Smith         col -= aij->cstart;
4346fd7127cSSatish Balay         ierr = MatSetValues(aij->A,1,&row,1,&col,&val,addv); CHKERRQ(ierr);
4353a40ed3dSBarry Smith       } else {
43655a1b374SBarry Smith         if (mat->was_assembled || mat->assembled) {
437905e6a2fSBarry Smith           if (!aij->colmap) {
438905e6a2fSBarry Smith             ierr = CreateColmap_MPIAIJ_Private(mat); CHKERRQ(ierr);
439905e6a2fSBarry Smith           }
440905e6a2fSBarry Smith           col = aij->colmap[col] - 1;
441ec8511deSBarry Smith           if (col < 0  && !((Mat_SeqAIJ*)(aij->A->data))->nonew) {
4422493cbb0SBarry Smith             ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
443227d817aSBarry Smith             col = (int) PetscReal(values[3*i+1]);
444d6dfbf8fSBarry Smith           }
4459e25ed09SBarry Smith         }
4466fd7127cSSatish Balay         ierr = MatSetValues(aij->B,1,&row,1,&col,&val,addv); CHKERRQ(ierr);
4471eb62cbbSBarry Smith       }
4481eb62cbbSBarry Smith     }
4491eb62cbbSBarry Smith     count--;
4501eb62cbbSBarry Smith   }
4510452661fSBarry Smith   PetscFree(aij->recv_waits); PetscFree(aij->rvalues);
4521eb62cbbSBarry Smith 
4531eb62cbbSBarry Smith   /* wait on sends */
4541eb62cbbSBarry Smith   if (aij->nsends) {
4550a198c4cSBarry Smith     send_status = (MPI_Status *) PetscMalloc(aij->nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
456ca161407SBarry Smith     ierr        = MPI_Waitall(aij->nsends,aij->send_waits,send_status);CHKERRQ(ierr);
4570452661fSBarry Smith     PetscFree(send_status);
4581eb62cbbSBarry Smith   }
459*b5bd3ad5SBarry Smith   if (aij->send_waits) PetscFree(aij->send_waits);
460*b5bd3ad5SBarry Smith   if (aij->svalues)    PetscFree(aij->svalues);
4611eb62cbbSBarry Smith 
46278b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr);
46378b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr);
4641eb62cbbSBarry Smith 
4652493cbb0SBarry Smith   /* determine if any processor has disassembled, if so we must
4662493cbb0SBarry Smith      also disassemble ourselfs, in order that we may reassemble. */
467ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);CHKERRQ(ierr);
468227d817aSBarry Smith   if (mat->was_assembled && !other_disassembled) {
4692493cbb0SBarry Smith     ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
4702493cbb0SBarry Smith   }
4712493cbb0SBarry Smith 
4726d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
47378b31e54SBarry Smith     ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr);
4745e42470aSBarry Smith   }
47578b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr);
47678b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr);
4775e42470aSBarry Smith 
4787a0afa10SBarry Smith   if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;}
4793a40ed3dSBarry Smith   PetscFunctionReturn(0);
4808a729477SBarry Smith }
4818a729477SBarry Smith 
4825615d1e5SSatish Balay #undef __FUNC__
4835615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIAIJ"
4848f6be9afSLois Curfman McInnes int MatZeroEntries_MPIAIJ(Mat A)
4851eb62cbbSBarry Smith {
48644a69424SLois Curfman McInnes   Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data;
487dbd7a890SLois Curfman McInnes   int        ierr;
4883a40ed3dSBarry Smith 
4893a40ed3dSBarry Smith   PetscFunctionBegin;
49078b31e54SBarry Smith   ierr = MatZeroEntries(l->A); CHKERRQ(ierr);
49178b31e54SBarry Smith   ierr = MatZeroEntries(l->B); CHKERRQ(ierr);
4923a40ed3dSBarry Smith   PetscFunctionReturn(0);
4931eb62cbbSBarry Smith }
4941eb62cbbSBarry Smith 
4955615d1e5SSatish Balay #undef __FUNC__
4965615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ"
4978f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag)
4981eb62cbbSBarry Smith {
49944a69424SLois Curfman McInnes   Mat_MPIAIJ     *l = (Mat_MPIAIJ *) A->data;
50017699dbbSLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
5016a5c57faSSatish Balay   int            *procs,*nprocs,j,found,idx,nsends,*work,row;
50217699dbbSLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
5035392566eSBarry Smith   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
5046a5c57faSSatish Balay   int            *lens,imdex,*lrows,*values,rstart=l->rstart;
505d6dfbf8fSBarry Smith   MPI_Comm       comm = A->comm;
5061eb62cbbSBarry Smith   MPI_Request    *send_waits,*recv_waits;
5071eb62cbbSBarry Smith   MPI_Status     recv_status,*send_status;
5081eb62cbbSBarry Smith   IS             istmp;
5096eb55b6aSBarry Smith   PetscTruth     localdiag;
5101eb62cbbSBarry Smith 
5113a40ed3dSBarry Smith   PetscFunctionBegin;
51277c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
51378b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5141eb62cbbSBarry Smith 
5151eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
5160452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
517cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
5180452661fSBarry Smith   owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
5191eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
5201eb62cbbSBarry Smith     idx = rows[i];
5211eb62cbbSBarry Smith     found = 0;
52217699dbbSLois Curfman McInnes     for ( j=0; j<size; j++ ) {
5231eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
5241eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
5251eb62cbbSBarry Smith       }
5261eb62cbbSBarry Smith     }
527a8c6a408SBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range");
5281eb62cbbSBarry Smith   }
52917699dbbSLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
5301eb62cbbSBarry Smith 
5311eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
5320452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
533ca161407SBarry Smith   ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
53417699dbbSLois Curfman McInnes   nrecvs = work[rank];
535ca161407SBarry Smith   ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
53617699dbbSLois Curfman McInnes   nmax = work[rank];
5370452661fSBarry Smith   PetscFree(work);
5381eb62cbbSBarry Smith 
5391eb62cbbSBarry Smith   /* post receives:   */
5403a40ed3dSBarry Smith   rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues);
541ca161407SBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits);
5421eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
543ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
5441eb62cbbSBarry Smith   }
5451eb62cbbSBarry Smith 
5461eb62cbbSBarry Smith   /* do sends:
5471eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
5481eb62cbbSBarry Smith          the ith processor
5491eb62cbbSBarry Smith   */
5500452661fSBarry Smith   svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
5513a40ed3dSBarry Smith   send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits);
5520452661fSBarry Smith   starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts);
5531eb62cbbSBarry Smith   starts[0] = 0;
55417699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
5551eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
5561eb62cbbSBarry Smith     svalues[starts[owner[i]]++] = rows[i];
5571eb62cbbSBarry Smith   }
5581eb62cbbSBarry Smith   ISRestoreIndices(is,&rows);
5591eb62cbbSBarry Smith 
5601eb62cbbSBarry Smith   starts[0] = 0;
56117699dbbSLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
5621eb62cbbSBarry Smith   count = 0;
56317699dbbSLois Curfman McInnes   for ( i=0; i<size; i++ ) {
5641eb62cbbSBarry Smith     if (procs[i]) {
565ca161407SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
5661eb62cbbSBarry Smith     }
5671eb62cbbSBarry Smith   }
5680452661fSBarry Smith   PetscFree(starts);
5691eb62cbbSBarry Smith 
57017699dbbSLois Curfman McInnes   base = owners[rank];
5711eb62cbbSBarry Smith 
5721eb62cbbSBarry Smith   /*  wait on receives */
5730452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
5741eb62cbbSBarry Smith   source = lens + nrecvs;
5751eb62cbbSBarry Smith   count  = nrecvs; slen = 0;
5761eb62cbbSBarry Smith   while (count) {
577ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
5781eb62cbbSBarry Smith     /* unpack receives into our local space */
579ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
580d6dfbf8fSBarry Smith     source[imdex]  = recv_status.MPI_SOURCE;
581d6dfbf8fSBarry Smith     lens[imdex]  = n;
5821eb62cbbSBarry Smith     slen += n;
5831eb62cbbSBarry Smith     count--;
5841eb62cbbSBarry Smith   }
5850452661fSBarry Smith   PetscFree(recv_waits);
5861eb62cbbSBarry Smith 
5871eb62cbbSBarry Smith   /* move the data into the send scatter */
5880452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
5891eb62cbbSBarry Smith   count = 0;
5901eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
5911eb62cbbSBarry Smith     values = rvalues + i*nmax;
5921eb62cbbSBarry Smith     for ( j=0; j<lens[i]; j++ ) {
5931eb62cbbSBarry Smith       lrows[count++] = values[j] - base;
5941eb62cbbSBarry Smith     }
5951eb62cbbSBarry Smith   }
5960452661fSBarry Smith   PetscFree(rvalues); PetscFree(lens);
5970452661fSBarry Smith   PetscFree(owner); PetscFree(nprocs);
5981eb62cbbSBarry Smith 
5991eb62cbbSBarry Smith   /* actually zap the local rows */
600029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
601464493b3SBarry Smith   PLogObjectParent(A,istmp);
6026a5c57faSSatish Balay 
6036eb55b6aSBarry Smith   /*
6046eb55b6aSBarry Smith         Zero the required rows. If the "diagonal block" of the matrix
6056eb55b6aSBarry Smith      is square and the user wishes to set the diagonal we use seperate
6066eb55b6aSBarry Smith      code so that MatSetValues() is not called for each diagonal allocating
6076eb55b6aSBarry Smith      new memory, thus calling lots of mallocs and slowing things down.
6086eb55b6aSBarry Smith 
6096eb55b6aSBarry Smith        Contributed by: Mathew Knepley
6106eb55b6aSBarry Smith   */
6116eb55b6aSBarry Smith   localdiag = PETSC_FALSE;
6126eb55b6aSBarry Smith   if (diag && (l->A->M == l->A->N)) {
6136eb55b6aSBarry Smith     localdiag = PETSC_TRUE;
6146eb55b6aSBarry Smith     ierr      = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
6156eb55b6aSBarry Smith   } else {
6166a5c57faSSatish Balay     ierr = MatZeroRows(l->A,istmp,0); CHKERRQ(ierr);
6176eb55b6aSBarry Smith   }
61878b31e54SBarry Smith   ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr);
61978b31e54SBarry Smith   ierr = ISDestroy(istmp); CHKERRQ(ierr);
6201eb62cbbSBarry Smith 
6216eb55b6aSBarry Smith   if (diag && (localdiag == PETSC_FALSE)) {
6226eb55b6aSBarry Smith     for ( i = 0; i < slen; i++ ) {
6236eb55b6aSBarry Smith       row = lrows[i] + rstart;
6246eb55b6aSBarry Smith       MatSetValues(A,1,&row,1,&row,diag,INSERT_VALUES);
6256eb55b6aSBarry Smith     }
6266eb55b6aSBarry Smith     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
6276eb55b6aSBarry Smith     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
6286eb55b6aSBarry Smith   }
6296eb55b6aSBarry Smith 
6306a5c57faSSatish Balay   if (diag) {
6316a5c57faSSatish Balay     for ( i = 0; i < slen; i++ ) {
6326a5c57faSSatish Balay       row = lrows[i] + rstart;
6336a5c57faSSatish Balay       MatSetValues(A,1,&row,1,&row,diag,INSERT_VALUES);
6346a5c57faSSatish Balay     }
6356a5c57faSSatish Balay     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
6366a5c57faSSatish Balay     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
6376a5c57faSSatish Balay   }
6386a5c57faSSatish Balay   PetscFree(lrows);
63972dacd9aSBarry Smith 
6401eb62cbbSBarry Smith   /* wait on sends */
6411eb62cbbSBarry Smith   if (nsends) {
642ca161407SBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status);
643ca161407SBarry Smith     ierr        = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
6440452661fSBarry Smith     PetscFree(send_status);
6451eb62cbbSBarry Smith   }
6460452661fSBarry Smith   PetscFree(send_waits); PetscFree(svalues);
6471eb62cbbSBarry Smith 
6483a40ed3dSBarry Smith   PetscFunctionReturn(0);
6491eb62cbbSBarry Smith }
6501eb62cbbSBarry Smith 
6515615d1e5SSatish Balay #undef __FUNC__
6525615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ"
6538f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy)
6541eb62cbbSBarry Smith {
655416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
656fbd6ef76SBarry Smith   int        ierr,nt;
657416022c9SBarry Smith 
6583a40ed3dSBarry Smith   PetscFunctionBegin;
659a2ce50c7SBarry Smith   ierr = VecGetLocalSize(xx,&nt);  CHKERRQ(ierr);
660fbd6ef76SBarry Smith   if (nt != a->n) {
661a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Incompatible partition of A and xx");
662fbd6ef76SBarry Smith   }
66343a90d84SBarry Smith   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr);
664f830108cSBarry Smith   ierr = (*a->A->ops->mult)(a->A,xx,yy); CHKERRQ(ierr);
66543a90d84SBarry Smith   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr);
666f830108cSBarry Smith   ierr = (*a->B->ops->multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr);
6673a40ed3dSBarry Smith   PetscFunctionReturn(0);
6681eb62cbbSBarry Smith }
6691eb62cbbSBarry Smith 
6705615d1e5SSatish Balay #undef __FUNC__
6715615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ"
6728f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
673da3a660dSBarry Smith {
674416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
675da3a660dSBarry Smith   int        ierr;
6763a40ed3dSBarry Smith 
6773a40ed3dSBarry Smith   PetscFunctionBegin;
67843a90d84SBarry Smith   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
679f830108cSBarry Smith   ierr = (*a->A->ops->multadd)(a->A,xx,yy,zz); CHKERRQ(ierr);
68043a90d84SBarry Smith   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
681f830108cSBarry Smith   ierr = (*a->B->ops->multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr);
6823a40ed3dSBarry Smith   PetscFunctionReturn(0);
683da3a660dSBarry Smith }
684da3a660dSBarry Smith 
6855615d1e5SSatish Balay #undef __FUNC__
6865615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ"
6878f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy)
688da3a660dSBarry Smith {
689416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
690da3a660dSBarry Smith   int        ierr;
691da3a660dSBarry Smith 
6923a40ed3dSBarry Smith   PetscFunctionBegin;
693da3a660dSBarry Smith   /* do nondiagonal part */
694f830108cSBarry Smith   ierr = (*a->B->ops->multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr);
695da3a660dSBarry Smith   /* send it on its way */
696537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
697da3a660dSBarry Smith   /* do local part */
698f830108cSBarry Smith   ierr = (*a->A->ops->multtrans)(a->A,xx,yy); CHKERRQ(ierr);
699da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
700da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
701da3a660dSBarry Smith   /* but is not perhaps always true. */
702537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
7033a40ed3dSBarry Smith   PetscFunctionReturn(0);
704da3a660dSBarry Smith }
705da3a660dSBarry Smith 
7065615d1e5SSatish Balay #undef __FUNC__
7075615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ"
7088f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
709da3a660dSBarry Smith {
710416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
711da3a660dSBarry Smith   int        ierr;
712da3a660dSBarry Smith 
7133a40ed3dSBarry Smith   PetscFunctionBegin;
714da3a660dSBarry Smith   /* do nondiagonal part */
715f830108cSBarry Smith   ierr = (*a->B->ops->multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr);
716da3a660dSBarry Smith   /* send it on its way */
717537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
718da3a660dSBarry Smith   /* do local part */
719f830108cSBarry Smith   ierr = (*a->A->ops->multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr);
720da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
721da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
722da3a660dSBarry Smith   /* but is not perhaps always true. */
7230a198c4cSBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
7243a40ed3dSBarry Smith   PetscFunctionReturn(0);
725da3a660dSBarry Smith }
726da3a660dSBarry Smith 
7271eb62cbbSBarry Smith /*
7281eb62cbbSBarry Smith   This only works correctly for square matrices where the subblock A->A is the
7291eb62cbbSBarry Smith    diagonal block
7301eb62cbbSBarry Smith */
7315615d1e5SSatish Balay #undef __FUNC__
7325615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ"
7338f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v)
7341eb62cbbSBarry Smith {
7353a40ed3dSBarry Smith   int        ierr;
736416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
7373a40ed3dSBarry Smith 
7383a40ed3dSBarry Smith   PetscFunctionBegin;
739a8c6a408SBarry Smith   if (a->M != a->N) SETERRQ(PETSC_ERR_SUP,0,"Supports only square matrix where A->A is diag block");
7405baf8537SBarry Smith   if (a->rstart != a->cstart || a->rend != a->cend) {
741a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"row partition must equal col partition");
7423a40ed3dSBarry Smith   }
7433a40ed3dSBarry Smith   ierr = MatGetDiagonal(a->A,v);CHKERRQ(ierr);
7443a40ed3dSBarry Smith   PetscFunctionReturn(0);
7451eb62cbbSBarry Smith }
7461eb62cbbSBarry Smith 
7475615d1e5SSatish Balay #undef __FUNC__
7485615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ"
7498f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A)
750052efed2SBarry Smith {
751052efed2SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
752052efed2SBarry Smith   int        ierr;
7533a40ed3dSBarry Smith 
7543a40ed3dSBarry Smith   PetscFunctionBegin;
755052efed2SBarry Smith   ierr = MatScale(aa,a->A); CHKERRQ(ierr);
756052efed2SBarry Smith   ierr = MatScale(aa,a->B); CHKERRQ(ierr);
7573a40ed3dSBarry Smith   PetscFunctionReturn(0);
758052efed2SBarry Smith }
759052efed2SBarry Smith 
7605615d1e5SSatish Balay #undef __FUNC__
761d4bb536fSBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ"
762e1311b90SBarry Smith int MatDestroy_MPIAIJ(Mat mat)
7631eb62cbbSBarry Smith {
76444a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
7651eb62cbbSBarry Smith   int        ierr;
76683e2fdc7SBarry Smith 
7673a40ed3dSBarry Smith   PetscFunctionBegin;
7683a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
769e1311b90SBarry Smith   PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",aij->M,aij->N);
770a5a9c739SBarry Smith #endif
77183e2fdc7SBarry Smith   ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr);
7720452661fSBarry Smith   PetscFree(aij->rowners);
77378b31e54SBarry Smith   ierr = MatDestroy(aij->A); CHKERRQ(ierr);
77478b31e54SBarry Smith   ierr = MatDestroy(aij->B); CHKERRQ(ierr);
7750452661fSBarry Smith   if (aij->colmap) PetscFree(aij->colmap);
7760452661fSBarry Smith   if (aij->garray) PetscFree(aij->garray);
7771eb62cbbSBarry Smith   if (aij->lvec)   VecDestroy(aij->lvec);
778a56f8943SBarry Smith   if (aij->Mvctx)  VecScatterDestroy(aij->Mvctx);
7797a0afa10SBarry Smith   if (aij->rowvalues) PetscFree(aij->rowvalues);
7800452661fSBarry Smith   PetscFree(aij);
781a5a9c739SBarry Smith   PLogObjectDestroy(mat);
7820452661fSBarry Smith   PetscHeaderDestroy(mat);
7833a40ed3dSBarry Smith   PetscFunctionReturn(0);
7841eb62cbbSBarry Smith }
785ee50ffe9SBarry Smith 
7865615d1e5SSatish Balay #undef __FUNC__
787d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary"
7888f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer)
7891eb62cbbSBarry Smith {
790416022c9SBarry Smith   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
791416022c9SBarry Smith   int         ierr;
792416022c9SBarry Smith 
7933a40ed3dSBarry Smith   PetscFunctionBegin;
79417699dbbSLois Curfman McInnes   if (aij->size == 1) {
795416022c9SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
796416022c9SBarry Smith   }
797a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported");
7983a40ed3dSBarry Smith   PetscFunctionReturn(0);
799416022c9SBarry Smith }
800416022c9SBarry Smith 
8015615d1e5SSatish Balay #undef __FUNC__
802d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworMatlab"
8038f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer)
804416022c9SBarry Smith {
80544a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
806dbb450caSBarry Smith   Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data;
807a56f8943SBarry Smith   int         ierr, format,shift = C->indexshift,rank;
808d636dbe3SBarry Smith   FILE        *fd;
80919bcc07fSBarry Smith   ViewerType  vtype;
810416022c9SBarry Smith 
8113a40ed3dSBarry Smith   PetscFunctionBegin;
81219bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
81319bcc07fSBarry Smith   if (vtype  == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) {
81490ace30eSBarry Smith     ierr = ViewerGetFormat(viewer,&format);
8150a198c4cSBarry Smith     if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
8164e220ebcSLois Curfman McInnes       MatInfo info;
8174e220ebcSLois Curfman McInnes       int     flg;
818a56f8943SBarry Smith       MPI_Comm_rank(mat->comm,&rank);
81990ace30eSBarry Smith       ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
8204e220ebcSLois Curfman McInnes       ierr = MatGetInfo(mat,MAT_LOCAL,&info);
82195e01e2fSLois Curfman McInnes       ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr);
82277c4ece6SBarry Smith       PetscSequentialPhaseBegin(mat->comm,1);
82395e01e2fSLois Curfman McInnes       if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n",
8244e220ebcSLois Curfman McInnes          rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
82595e01e2fSLois Curfman McInnes       else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n",
8264e220ebcSLois Curfman McInnes          rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
8274e220ebcSLois Curfman McInnes       ierr = MatGetInfo(aij->A,MAT_LOCAL,&info);
8284e220ebcSLois Curfman McInnes       fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used);
8294e220ebcSLois Curfman McInnes       ierr = MatGetInfo(aij->B,MAT_LOCAL,&info);
8304e220ebcSLois Curfman McInnes       fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used);
831a56f8943SBarry Smith       fflush(fd);
83277c4ece6SBarry Smith       PetscSequentialPhaseEnd(mat->comm,1);
833a40aa06bSLois Curfman McInnes       ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr);
8343a40ed3dSBarry Smith       PetscFunctionReturn(0);
8353a40ed3dSBarry Smith     } else if (format == VIEWER_FORMAT_ASCII_INFO) {
8363a40ed3dSBarry Smith       PetscFunctionReturn(0);
83708480c60SBarry Smith     }
838416022c9SBarry Smith   }
839416022c9SBarry Smith 
84019bcc07fSBarry Smith   if (vtype == DRAW_VIEWER) {
84119bcc07fSBarry Smith     Draw       draw;
84219bcc07fSBarry Smith     PetscTruth isnull;
84319bcc07fSBarry Smith     ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
8443a40ed3dSBarry Smith     ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
84519bcc07fSBarry Smith   }
84619bcc07fSBarry Smith 
84719bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
84890ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
84977c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
850d13ab20cSBarry Smith     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n",
85117699dbbSLois Curfman McInnes            aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart,
8521eb62cbbSBarry Smith            aij->cend);
85378b31e54SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
85478b31e54SBarry Smith     ierr = MatView(aij->B,viewer); CHKERRQ(ierr);
855d13ab20cSBarry Smith     fflush(fd);
85677c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
8573a40ed3dSBarry Smith   } else {
858a56f8943SBarry Smith     int size = aij->size;
859a56f8943SBarry Smith     rank = aij->rank;
86017699dbbSLois Curfman McInnes     if (size == 1) {
86178b31e54SBarry Smith       ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
8623a40ed3dSBarry Smith     } else {
86395373324SBarry Smith       /* assemble the entire matrix onto first processor. */
86495373324SBarry Smith       Mat         A;
865ec8511deSBarry Smith       Mat_SeqAIJ *Aloc;
8662eb8c8abSBarry Smith       int         M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct;
86795373324SBarry Smith       Scalar      *a;
8682ee70a88SLois Curfman McInnes 
86917699dbbSLois Curfman McInnes       if (!rank) {
87055843e3eSBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr);
8713a40ed3dSBarry Smith       } else {
87255843e3eSBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);CHKERRQ(ierr);
87395373324SBarry Smith       }
874464493b3SBarry Smith       PLogObjectParent(mat,A);
875416022c9SBarry Smith 
87695373324SBarry Smith       /* copy over the A part */
877ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->A->data;
8782ee70a88SLois Curfman McInnes       m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
87995373324SBarry Smith       row = aij->rstart;
880dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;}
88195373324SBarry Smith       for ( i=0; i<m; i++ ) {
882416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr);
88395373324SBarry Smith         row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
88495373324SBarry Smith       }
8852ee70a88SLois Curfman McInnes       aj = Aloc->j;
886dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;}
88795373324SBarry Smith 
88895373324SBarry Smith       /* copy over the B part */
889ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->B->data;
8902ee70a88SLois Curfman McInnes       m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
89195373324SBarry Smith       row = aij->rstart;
8920452661fSBarry Smith       ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
893dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];}
89495373324SBarry Smith       for ( i=0; i<m; i++ ) {
895416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr);
89695373324SBarry Smith         row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
89795373324SBarry Smith       }
8980452661fSBarry Smith       PetscFree(ct);
8996d4a8577SBarry Smith       ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9006d4a8577SBarry Smith       ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
90155843e3eSBarry Smith       /*
90255843e3eSBarry Smith          Everyone has to call to draw the matrix since the graphics waits are
90355843e3eSBarry Smith          synchronized across all processors that share the Draw object
90455843e3eSBarry Smith       */
90555843e3eSBarry Smith       if (!rank || vtype == DRAW_VIEWER) {
90678b31e54SBarry Smith         ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr);
90795373324SBarry Smith       }
90878b31e54SBarry Smith       ierr = MatDestroy(A); CHKERRQ(ierr);
90995373324SBarry Smith     }
91095373324SBarry Smith   }
9113a40ed3dSBarry Smith   PetscFunctionReturn(0);
9121eb62cbbSBarry Smith }
9131eb62cbbSBarry Smith 
9145615d1e5SSatish Balay #undef __FUNC__
915d4bb536fSBarry Smith #define __FUNC__ "MatView_MPIAIJ"
916e1311b90SBarry Smith int MatView_MPIAIJ(Mat mat,Viewer viewer)
917416022c9SBarry Smith {
918416022c9SBarry Smith   int         ierr;
91919bcc07fSBarry Smith   ViewerType  vtype;
920416022c9SBarry Smith 
9213a40ed3dSBarry Smith   PetscFunctionBegin;
92219bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
92319bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER ||
92419bcc07fSBarry Smith       vtype == DRAW_VIEWER       || vtype == MATLAB_VIEWER) {
925d7e8b826SBarry Smith     ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr);
9265cd90555SBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
9273a40ed3dSBarry Smith     ierr = MatView_MPIAIJ_Binary(mat,viewer);CHKERRQ(ierr);
9285cd90555SBarry Smith   } else {
9295cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
930416022c9SBarry Smith   }
9313a40ed3dSBarry Smith   PetscFunctionReturn(0);
932416022c9SBarry Smith }
933416022c9SBarry Smith 
9341eb62cbbSBarry Smith /*
9351eb62cbbSBarry Smith     This has to provide several versions.
9361eb62cbbSBarry Smith 
9371eb62cbbSBarry Smith      2) a) use only local smoothing updating outer values only once.
9381eb62cbbSBarry Smith         b) local smoothing updating outer values each inner iteration
939d6dfbf8fSBarry Smith      3) color updating out values betwen colors.
9401eb62cbbSBarry Smith */
9415615d1e5SSatish Balay #undef __FUNC__
9425615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ"
9438f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag,
944dbb450caSBarry Smith                            double fshift,int its,Vec xx)
9458a729477SBarry Smith {
94644a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
947d6dfbf8fSBarry Smith   Mat        AA = mat->A, BB = mat->B;
948ec8511deSBarry Smith   Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data;
949c16cb8f2SBarry Smith   Scalar     *b,*x,*xs,*ls,d,*v,sum;
9506abc6512SBarry Smith   int        ierr,*idx, *diag;
951416022c9SBarry Smith   int        n = mat->n, m = mat->m, i,shift = A->indexshift;
9528a729477SBarry Smith 
9533a40ed3dSBarry Smith   PetscFunctionBegin;
954d6dfbf8fSBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls);
955dbb450caSBarry Smith   xs = x + shift; /* shift by one for index start of 1 */
956dbb450caSBarry Smith   ls = ls + shift;
95783e2fdc7SBarry Smith   if (!A->diag) {ierr = MatMarkDiag_SeqAIJ(AA); CHKERRQ(ierr);}
958d6dfbf8fSBarry Smith   diag = A->diag;
959c16cb8f2SBarry Smith   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
960da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
961f830108cSBarry Smith       ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr);
9623a40ed3dSBarry Smith       PetscFunctionReturn(0);
963da3a660dSBarry Smith     }
9643a40ed3dSBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr);
9653a40ed3dSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr);
966d6dfbf8fSBarry Smith     while (its--) {
967d6dfbf8fSBarry Smith       /* go down through the rows */
968d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
969d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
970dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
971dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
972d6dfbf8fSBarry Smith         sum  = b[i];
973d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
974dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
975d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
976dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
977dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
978d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
97955a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
980d6dfbf8fSBarry Smith       }
981d6dfbf8fSBarry Smith       /* come up through the rows */
982d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
983d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
984dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
985dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
986d6dfbf8fSBarry Smith         sum  = b[i];
987d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
988dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
989d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
990dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
991dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
992d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
99355a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
994d6dfbf8fSBarry Smith       }
995d6dfbf8fSBarry Smith     }
9963a40ed3dSBarry Smith   } else if (flag & SOR_LOCAL_FORWARD_SWEEP){
997da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
998f830108cSBarry Smith       ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr);
9993a40ed3dSBarry Smith       PetscFunctionReturn(0);
1000da3a660dSBarry Smith     }
10013a40ed3dSBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr);
10023a40ed3dSBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);CHKERRQ(ierr);
1003d6dfbf8fSBarry Smith     while (its--) {
1004d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
1005d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
1006dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
1007dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
1008d6dfbf8fSBarry Smith         sum  = b[i];
1009d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1010dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
1011d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
1012dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
1013dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
1014d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
101555a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
1016d6dfbf8fSBarry Smith       }
1017d6dfbf8fSBarry Smith     }
10183a40ed3dSBarry Smith   } else if (flag & SOR_LOCAL_BACKWARD_SWEEP){
1019da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
1020f830108cSBarry Smith       ierr = (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,its,xx);CHKERRQ(ierr);
10213a40ed3dSBarry Smith       PetscFunctionReturn(0);
1022da3a660dSBarry Smith     }
102343a90d84SBarry Smith     ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,
102478b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
102543a90d84SBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,
102678b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
1027d6dfbf8fSBarry Smith     while (its--) {
1028d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
1029d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
1030dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
1031dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
1032d6dfbf8fSBarry Smith         sum  = b[i];
1033d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1034dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
1035d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
1036dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
1037dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
1038d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
103955a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
1040d6dfbf8fSBarry Smith       }
1041d6dfbf8fSBarry Smith     }
10423a40ed3dSBarry Smith   } else {
1043a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Parallel SOR not supported");
1044c16cb8f2SBarry Smith   }
10453a40ed3dSBarry Smith   PetscFunctionReturn(0);
10468a729477SBarry Smith }
1047a66be287SLois Curfman McInnes 
10485615d1e5SSatish Balay #undef __FUNC__
1049d4bb536fSBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ"
10508f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1051a66be287SLois Curfman McInnes {
1052a66be287SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1053a66be287SLois Curfman McInnes   Mat        A = mat->A, B = mat->B;
10544e220ebcSLois Curfman McInnes   int        ierr;
10554e220ebcSLois Curfman McInnes   double     isend[5], irecv[5];
1056a66be287SLois Curfman McInnes 
10573a40ed3dSBarry Smith   PetscFunctionBegin;
10584e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10594e220ebcSLois Curfman McInnes   ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr);
10604e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
10614e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
10624e220ebcSLois Curfman McInnes   ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr);
10634e220ebcSLois Curfman McInnes   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
10644e220ebcSLois Curfman McInnes   isend[3] += info->memory;  isend[4] += info->mallocs;
1065a66be287SLois Curfman McInnes   if (flag == MAT_LOCAL) {
10664e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
10674e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
10684e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
10694e220ebcSLois Curfman McInnes     info->memory       = isend[3];
10704e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
1071a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
1072ca161407SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm);CHKERRQ(ierr);
10734e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
10744e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
10754e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
10764e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
10774e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
1078a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
1079ca161407SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm);CHKERRQ(ierr);
10804e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
10814e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
10824e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
10834e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
10844e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
1085a66be287SLois Curfman McInnes   }
10864e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
10874e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
10884e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
10898d700155SBarry Smith   info->rows_global       = (double)mat->M;
10908d700155SBarry Smith   info->columns_global    = (double)mat->N;
10918d700155SBarry Smith   info->rows_local        = (double)mat->m;
10928d700155SBarry Smith   info->columns_local     = (double)mat->N;
10934e220ebcSLois Curfman McInnes 
10943a40ed3dSBarry Smith   PetscFunctionReturn(0);
1095a66be287SLois Curfman McInnes }
1096a66be287SLois Curfman McInnes 
10975615d1e5SSatish Balay #undef __FUNC__
1098d4bb536fSBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ"
10998f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op)
1100c74985f6SBarry Smith {
1101c0bbcb79SLois Curfman McInnes   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1102c74985f6SBarry Smith 
11033a40ed3dSBarry Smith   PetscFunctionBegin;
11046d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
11056d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
11066da5968aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED ||
1107c2653b3dSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
110896854ed6SLois Curfman McInnes       op == MAT_NEW_NONZERO_ALLOCATION_ERROR ||
1109c2653b3dSLois Curfman McInnes       op == MAT_NEW_NONZERO_LOCATION_ERROR) {
1110b1fbbac0SLois Curfman McInnes         MatSetOption(a->A,op);
1111b1fbbac0SLois Curfman McInnes         MatSetOption(a->B,op);
1112b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
1113aeafbbfcSLois Curfman McInnes     a->roworiented = 1;
1114c0bbcb79SLois Curfman McInnes     MatSetOption(a->A,op);
1115c0bbcb79SLois Curfman McInnes     MatSetOption(a->B,op);
1116b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
11176da5968aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
11186d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
11196d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
11206d4a8577SBarry Smith              op == MAT_YES_NEW_DIAGONALS)
1121981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIAIJ:Option ignored\n");
11226d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED) {
11234b0e389bSBarry Smith     a->roworiented = 0;
11244b0e389bSBarry Smith     MatSetOption(a->A,op);
11254b0e389bSBarry Smith     MatSetOption(a->B,op);
11262b362799SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
112790f02eecSBarry Smith     a->donotstash = 1;
11283a40ed3dSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS){
11293a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
11303a40ed3dSBarry Smith   } else {
11313a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
11323a40ed3dSBarry Smith   }
11333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1134c74985f6SBarry Smith }
1135c74985f6SBarry Smith 
11365615d1e5SSatish Balay #undef __FUNC__
1137d4bb536fSBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ"
11388f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n)
1139c74985f6SBarry Smith {
114044a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
11413a40ed3dSBarry Smith 
11423a40ed3dSBarry Smith   PetscFunctionBegin;
11430752156aSBarry Smith   if (m) *m = mat->M;
11440752156aSBarry Smith   if (n) *n = mat->N;
11453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1146c74985f6SBarry Smith }
1147c74985f6SBarry Smith 
11485615d1e5SSatish Balay #undef __FUNC__
1149d4bb536fSBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ"
11508f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n)
1151c74985f6SBarry Smith {
115244a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
11533a40ed3dSBarry Smith 
11543a40ed3dSBarry Smith   PetscFunctionBegin;
11550752156aSBarry Smith   if (m) *m = mat->m;
1156f830108cSBarry Smith   if (n) *n = mat->n;
11573a40ed3dSBarry Smith   PetscFunctionReturn(0);
1158c74985f6SBarry Smith }
1159c74985f6SBarry Smith 
11605615d1e5SSatish Balay #undef __FUNC__
1161d4bb536fSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ"
11628f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n)
1163c74985f6SBarry Smith {
116444a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
11653a40ed3dSBarry Smith 
11663a40ed3dSBarry Smith   PetscFunctionBegin;
1167c74985f6SBarry Smith   *m = mat->rstart; *n = mat->rend;
11683a40ed3dSBarry Smith   PetscFunctionReturn(0);
1169c74985f6SBarry Smith }
1170c74985f6SBarry Smith 
11716d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**);
11726d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**);
11736d84be18SBarry Smith 
11745615d1e5SSatish Balay #undef __FUNC__
11755615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ"
11766d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v)
117739e00950SLois Curfman McInnes {
1178154123eaSLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
117970f0671dSBarry Smith   Scalar     *vworkA, *vworkB, **pvA, **pvB,*v_p;
1180154123eaSLois Curfman McInnes   int        i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart;
1181154123eaSLois Curfman McInnes   int        nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend;
118270f0671dSBarry Smith   int        *cmap, *idx_p;
118339e00950SLois Curfman McInnes 
11843a40ed3dSBarry Smith   PetscFunctionBegin;
1185a8c6a408SBarry Smith   if (mat->getrowactive == PETSC_TRUE) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Already active");
11867a0afa10SBarry Smith   mat->getrowactive = PETSC_TRUE;
11877a0afa10SBarry Smith 
118870f0671dSBarry Smith   if (!mat->rowvalues && (idx || v)) {
11897a0afa10SBarry Smith     /*
11907a0afa10SBarry Smith         allocate enough space to hold information from the longest row.
11917a0afa10SBarry Smith     */
11927a0afa10SBarry Smith     Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data;
1193c16cb8f2SBarry Smith     int     max = 1,m = mat->m,tmp;
1194c16cb8f2SBarry Smith     for ( i=0; i<m; i++ ) {
11957a0afa10SBarry Smith       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i];
11967a0afa10SBarry Smith       if (max < tmp) { max = tmp; }
11977a0afa10SBarry Smith     }
11987a0afa10SBarry Smith     mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar)));
11997a0afa10SBarry Smith     CHKPTRQ(mat->rowvalues);
12007a0afa10SBarry Smith     mat->rowindices = (int *) (mat->rowvalues + max);
12017a0afa10SBarry Smith   }
12027a0afa10SBarry Smith 
1203a8c6a408SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Only local rows")
1204abc0e9e4SLois Curfman McInnes   lrow = row - rstart;
120539e00950SLois Curfman McInnes 
1206154123eaSLois Curfman McInnes   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1207154123eaSLois Curfman McInnes   if (!v)   {pvA = 0; pvB = 0;}
1208154123eaSLois Curfman McInnes   if (!idx) {pcA = 0; if (!v) pcB = 0;}
1209f830108cSBarry Smith   ierr = (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
1210f830108cSBarry Smith   ierr = (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
1211154123eaSLois Curfman McInnes   nztot = nzA + nzB;
1212154123eaSLois Curfman McInnes 
121370f0671dSBarry Smith   cmap  = mat->garray;
1214154123eaSLois Curfman McInnes   if (v  || idx) {
1215154123eaSLois Curfman McInnes     if (nztot) {
1216154123eaSLois Curfman McInnes       /* Sort by increasing column numbers, assuming A and B already sorted */
121770f0671dSBarry Smith       int imark = -1;
1218154123eaSLois Curfman McInnes       if (v) {
121970f0671dSBarry Smith         *v = v_p = mat->rowvalues;
122039e00950SLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
122170f0671dSBarry Smith           if (cmap[cworkB[i]] < cstart)   v_p[i] = vworkB[i];
1222154123eaSLois Curfman McInnes           else break;
1223154123eaSLois Curfman McInnes         }
1224154123eaSLois Curfman McInnes         imark = i;
122570f0671dSBarry Smith         for ( i=0; i<nzA; i++ )     v_p[imark+i] = vworkA[i];
122670f0671dSBarry Smith         for ( i=imark; i<nzB; i++ ) v_p[nzA+i]   = vworkB[i];
1227154123eaSLois Curfman McInnes       }
1228154123eaSLois Curfman McInnes       if (idx) {
122970f0671dSBarry Smith         *idx = idx_p = mat->rowindices;
123070f0671dSBarry Smith         if (imark > -1) {
123170f0671dSBarry Smith           for ( i=0; i<imark; i++ ) {
123270f0671dSBarry Smith             idx_p[i] = cmap[cworkB[i]];
123370f0671dSBarry Smith           }
123470f0671dSBarry Smith         } else {
1235154123eaSLois Curfman McInnes           for ( i=0; i<nzB; i++ ) {
123670f0671dSBarry Smith             if (cmap[cworkB[i]] < cstart)   idx_p[i] = cmap[cworkB[i]];
1237154123eaSLois Curfman McInnes             else break;
1238154123eaSLois Curfman McInnes           }
1239154123eaSLois Curfman McInnes           imark = i;
124070f0671dSBarry Smith         }
124170f0671dSBarry Smith         for ( i=0; i<nzA; i++ )     idx_p[imark+i] = cstart + cworkA[i];
124270f0671dSBarry Smith         for ( i=imark; i<nzB; i++ ) idx_p[nzA+i]   = cmap[cworkB[i]];
124339e00950SLois Curfman McInnes       }
124439e00950SLois Curfman McInnes     }
12451ca473b0SSatish Balay     else {
12461ca473b0SSatish Balay       if (idx) *idx = 0;
12471ca473b0SSatish Balay       if (v)   *v   = 0;
12481ca473b0SSatish Balay     }
1249154123eaSLois Curfman McInnes   }
125039e00950SLois Curfman McInnes   *nz = nztot;
1251f830108cSBarry Smith   ierr = (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
1252f830108cSBarry Smith   ierr = (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
12533a40ed3dSBarry Smith   PetscFunctionReturn(0);
125439e00950SLois Curfman McInnes }
125539e00950SLois Curfman McInnes 
12565615d1e5SSatish Balay #undef __FUNC__
1257d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ"
12586d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v)
125939e00950SLois Curfman McInnes {
12607a0afa10SBarry Smith   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
12613a40ed3dSBarry Smith 
12623a40ed3dSBarry Smith   PetscFunctionBegin;
12637a0afa10SBarry Smith   if (aij->getrowactive == PETSC_FALSE) {
1264a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"MatGetRow not called");
12657a0afa10SBarry Smith   }
12667a0afa10SBarry Smith   aij->getrowactive = PETSC_FALSE;
12673a40ed3dSBarry Smith   PetscFunctionReturn(0);
126839e00950SLois Curfman McInnes }
126939e00950SLois Curfman McInnes 
12705615d1e5SSatish Balay #undef __FUNC__
12715615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ"
12728f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm)
1273855ac2c5SLois Curfman McInnes {
1274855ac2c5SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1275ec8511deSBarry Smith   Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data;
1276416022c9SBarry Smith   int        ierr, i, j, cstart = aij->cstart,shift = amat->indexshift;
1277416022c9SBarry Smith   double     sum = 0.0;
127804ca555eSLois Curfman McInnes   Scalar     *v;
127904ca555eSLois Curfman McInnes 
12803a40ed3dSBarry Smith   PetscFunctionBegin;
128117699dbbSLois Curfman McInnes   if (aij->size == 1) {
128214183eadSLois Curfman McInnes     ierr =  MatNorm(aij->A,type,norm); CHKERRQ(ierr);
128337fa93a5SLois Curfman McInnes   } else {
128404ca555eSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
128504ca555eSLois Curfman McInnes       v = amat->a;
128604ca555eSLois Curfman McInnes       for (i=0; i<amat->nz; i++ ) {
12873a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
128804ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
128904ca555eSLois Curfman McInnes #else
129004ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
129104ca555eSLois Curfman McInnes #endif
129204ca555eSLois Curfman McInnes       }
129304ca555eSLois Curfman McInnes       v = bmat->a;
129404ca555eSLois Curfman McInnes       for (i=0; i<bmat->nz; i++ ) {
12953a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
129604ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
129704ca555eSLois Curfman McInnes #else
129804ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
129904ca555eSLois Curfman McInnes #endif
130004ca555eSLois Curfman McInnes       }
1301ca161407SBarry Smith       ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr);
130204ca555eSLois Curfman McInnes       *norm = sqrt(*norm);
13033a40ed3dSBarry Smith     } else if (type == NORM_1) { /* max column norm */
130404ca555eSLois Curfman McInnes       double *tmp, *tmp2;
130504ca555eSLois Curfman McInnes       int    *jj, *garray = aij->garray;
1306758f045eSSatish Balay       tmp  = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp);
1307758f045eSSatish Balay       tmp2 = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp2);
1308cddf8d76SBarry Smith       PetscMemzero(tmp,aij->N*sizeof(double));
130904ca555eSLois Curfman McInnes       *norm = 0.0;
131004ca555eSLois Curfman McInnes       v = amat->a; jj = amat->j;
131104ca555eSLois Curfman McInnes       for ( j=0; j<amat->nz; j++ ) {
1312579c6b6fSBarry Smith         tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v);  v++;
131304ca555eSLois Curfman McInnes       }
131404ca555eSLois Curfman McInnes       v = bmat->a; jj = bmat->j;
131504ca555eSLois Curfman McInnes       for ( j=0; j<bmat->nz; j++ ) {
1316579c6b6fSBarry Smith         tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++;
131704ca555eSLois Curfman McInnes       }
1318ca161407SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);CHKERRQ(ierr);
131904ca555eSLois Curfman McInnes       for ( j=0; j<aij->N; j++ ) {
132004ca555eSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
132104ca555eSLois Curfman McInnes       }
13220452661fSBarry Smith       PetscFree(tmp); PetscFree(tmp2);
13233a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
1324515d9167SLois Curfman McInnes       double ntemp = 0.0;
132504ca555eSLois Curfman McInnes       for ( j=0; j<amat->m; j++ ) {
1326dbb450caSBarry Smith         v = amat->a + amat->i[j] + shift;
132704ca555eSLois Curfman McInnes         sum = 0.0;
132804ca555eSLois Curfman McInnes         for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) {
1329cddf8d76SBarry Smith           sum += PetscAbsScalar(*v); v++;
133004ca555eSLois Curfman McInnes         }
1331dbb450caSBarry Smith         v = bmat->a + bmat->i[j] + shift;
133204ca555eSLois Curfman McInnes         for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) {
1333cddf8d76SBarry Smith           sum += PetscAbsScalar(*v); v++;
133404ca555eSLois Curfman McInnes         }
1335515d9167SLois Curfman McInnes         if (sum > ntemp) ntemp = sum;
133604ca555eSLois Curfman McInnes       }
1337ca161407SBarry Smith       ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);CHKERRQ(ierr);
1338ca161407SBarry Smith     } else {
1339a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
134004ca555eSLois Curfman McInnes     }
134137fa93a5SLois Curfman McInnes   }
13423a40ed3dSBarry Smith   PetscFunctionReturn(0);
1343855ac2c5SLois Curfman McInnes }
1344855ac2c5SLois Curfman McInnes 
13455615d1e5SSatish Balay #undef __FUNC__
13465615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ"
13478f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout)
1348b7c46309SBarry Smith {
1349b7c46309SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1350dbb450caSBarry Smith   Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data;
1351416022c9SBarry Smith   int        ierr,shift = Aloc->indexshift;
1352b7c46309SBarry Smith   int        M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct;
13533a40ed3dSBarry Smith   Mat        B;
1354b7c46309SBarry Smith   Scalar     *array;
1355b7c46309SBarry Smith 
13563a40ed3dSBarry Smith   PetscFunctionBegin;
1357d4bb536fSBarry Smith   if (matout == PETSC_NULL && M != N) {
1358a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place");
1359d4bb536fSBarry Smith   }
1360d4bb536fSBarry Smith 
1361d4bb536fSBarry Smith   ierr = MatCreateMPIAIJ(A->comm,a->n,a->m,N,M,0,PETSC_NULL,0,PETSC_NULL,&B);CHKERRQ(ierr);
1362b7c46309SBarry Smith 
1363b7c46309SBarry Smith   /* copy over the A part */
1364ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->A->data;
1365b7c46309SBarry Smith   m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1366b7c46309SBarry Smith   row = a->rstart;
1367dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;}
1368b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1369416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1370b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
1371b7c46309SBarry Smith   }
1372b7c46309SBarry Smith   aj = Aloc->j;
13734af08d9eSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;}
1374b7c46309SBarry Smith 
1375b7c46309SBarry Smith   /* copy over the B part */
1376ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->B->data;
1377b7c46309SBarry Smith   m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1378b7c46309SBarry Smith   row = a->rstart;
13790452661fSBarry Smith   ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols);
1380dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];}
1381b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1382416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1383b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
1384b7c46309SBarry Smith   }
13850452661fSBarry Smith   PetscFree(ct);
13866d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13876d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13883638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
13890de55854SLois Curfman McInnes     *matout = B;
13900de55854SLois Curfman McInnes   } else {
1391f830108cSBarry Smith     PetscOps       *Abops;
1392f830108cSBarry Smith     struct _MatOps *Aops;
1393f830108cSBarry Smith 
13940de55854SLois Curfman McInnes     /* This isn't really an in-place transpose .... but free data structures from a */
13950452661fSBarry Smith     PetscFree(a->rowners);
13960de55854SLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
13970de55854SLois Curfman McInnes     ierr = MatDestroy(a->B); CHKERRQ(ierr);
13980452661fSBarry Smith     if (a->colmap) PetscFree(a->colmap);
13990452661fSBarry Smith     if (a->garray) PetscFree(a->garray);
14000de55854SLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
1401a56f8943SBarry Smith     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
14020452661fSBarry Smith     PetscFree(a);
1403f830108cSBarry Smith 
1404f830108cSBarry Smith     /*
1405f830108cSBarry Smith        This is horrible, horrible code. We need to keep the
1406f830108cSBarry Smith       A pointers for the bops and ops but copy everything
1407f830108cSBarry Smith       else from C.
1408f830108cSBarry Smith     */
1409f830108cSBarry Smith     Abops = A->bops;
1410f830108cSBarry Smith     Aops  = A->ops;
1411f09e8eb9SSatish Balay     PetscMemcpy(A,B,sizeof(struct _p_Mat));
1412f830108cSBarry Smith     A->bops = Abops;
1413f830108cSBarry Smith     A->ops  = Aops;
14140452661fSBarry Smith     PetscHeaderDestroy(B);
14150de55854SLois Curfman McInnes   }
14163a40ed3dSBarry Smith   PetscFunctionReturn(0);
1417b7c46309SBarry Smith }
1418b7c46309SBarry Smith 
14195615d1e5SSatish Balay #undef __FUNC__
14205615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ"
14214b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr)
1422a008b906SSatish Balay {
14234b967eb1SSatish Balay   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
14244b967eb1SSatish Balay   Mat a = aij->A, b = aij->B;
1425a008b906SSatish Balay   int ierr,s1,s2,s3;
1426a008b906SSatish Balay 
14273a40ed3dSBarry Smith   PetscFunctionBegin;
14284b967eb1SSatish Balay   ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr);
14294b967eb1SSatish Balay   if (rr) {
1430e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&s1);CHKERRQ(ierr);
1431a8c6a408SBarry Smith     if (s1!=s3) SETERRQ(PETSC_ERR_ARG_SIZ,0,"right vector non-conforming local size");
14324b967eb1SSatish Balay     /* Overlap communication with computation. */
143343a90d84SBarry Smith     ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr);
1434a008b906SSatish Balay   }
14354b967eb1SSatish Balay   if (ll) {
1436e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&s1);CHKERRQ(ierr);
1437a8c6a408SBarry Smith     if (s1!=s2) SETERRQ(PETSC_ERR_ARG_SIZ,0,"left vector non-conforming local size");
1438f830108cSBarry Smith     ierr = (*b->ops->diagonalscale)(b,ll,0); CHKERRQ(ierr);
14394b967eb1SSatish Balay   }
14404b967eb1SSatish Balay   /* scale  the diagonal block */
1441f830108cSBarry Smith   ierr = (*a->ops->diagonalscale)(a,ll,rr); CHKERRQ(ierr);
14424b967eb1SSatish Balay 
14434b967eb1SSatish Balay   if (rr) {
14444b967eb1SSatish Balay     /* Do a scatter end and then right scale the off-diagonal block */
144543a90d84SBarry Smith     ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr);
1446f830108cSBarry Smith     ierr = (*b->ops->diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr);
14474b967eb1SSatish Balay   }
14484b967eb1SSatish Balay 
14493a40ed3dSBarry Smith   PetscFunctionReturn(0);
1450a008b906SSatish Balay }
1451a008b906SSatish Balay 
1452a008b906SSatish Balay 
1453682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat);
14545615d1e5SSatish Balay #undef __FUNC__
1455d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ"
14568f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A)
1457682d7d0cSBarry Smith {
1458682d7d0cSBarry Smith   Mat_MPIAIJ *a   = (Mat_MPIAIJ*) A->data;
14593a40ed3dSBarry Smith   int        ierr;
1460682d7d0cSBarry Smith 
14613a40ed3dSBarry Smith   PetscFunctionBegin;
14623a40ed3dSBarry Smith   if (!a->rank) {
14633a40ed3dSBarry Smith     ierr = MatPrintHelp_SeqAIJ(a->A);CHKERRQ(ierr);
14643a40ed3dSBarry Smith   }
14653a40ed3dSBarry Smith   PetscFunctionReturn(0);
1466682d7d0cSBarry Smith }
1467682d7d0cSBarry Smith 
14685615d1e5SSatish Balay #undef __FUNC__
1469d4bb536fSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ"
14708f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs)
14715a838052SSatish Balay {
14723a40ed3dSBarry Smith   PetscFunctionBegin;
14735a838052SSatish Balay   *bs = 1;
14743a40ed3dSBarry Smith   PetscFunctionReturn(0);
14755a838052SSatish Balay }
14765615d1e5SSatish Balay #undef __FUNC__
1477d4bb536fSBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ"
14788f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A)
1479bb5a7306SBarry Smith {
1480bb5a7306SBarry Smith   Mat_MPIAIJ *a   = (Mat_MPIAIJ*) A->data;
1481bb5a7306SBarry Smith   int        ierr;
14823a40ed3dSBarry Smith 
14833a40ed3dSBarry Smith   PetscFunctionBegin;
1484bb5a7306SBarry Smith   ierr = MatSetUnfactored(a->A); CHKERRQ(ierr);
14853a40ed3dSBarry Smith   PetscFunctionReturn(0);
1486bb5a7306SBarry Smith }
1487bb5a7306SBarry Smith 
1488d4bb536fSBarry Smith #undef __FUNC__
1489d4bb536fSBarry Smith #define __FUNC__ "MatEqual_MPIAIJ"
1490d4bb536fSBarry Smith int MatEqual_MPIAIJ(Mat A, Mat B, PetscTruth *flag)
1491d4bb536fSBarry Smith {
1492d4bb536fSBarry Smith   Mat_MPIAIJ *matB = (Mat_MPIAIJ *) B->data,*matA = (Mat_MPIAIJ *) A->data;
1493d4bb536fSBarry Smith   Mat        a, b, c, d;
1494d4bb536fSBarry Smith   PetscTruth flg;
1495d4bb536fSBarry Smith   int        ierr;
1496d4bb536fSBarry Smith 
14973a40ed3dSBarry Smith   PetscFunctionBegin;
1498a8c6a408SBarry Smith   if (B->type != MATMPIAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
1499d4bb536fSBarry Smith   a = matA->A; b = matA->B;
1500d4bb536fSBarry Smith   c = matB->A; d = matB->B;
1501d4bb536fSBarry Smith 
1502d4bb536fSBarry Smith   ierr = MatEqual(a, c, &flg); CHKERRQ(ierr);
1503d4bb536fSBarry Smith   if (flg == PETSC_TRUE) {
1504d4bb536fSBarry Smith     ierr = MatEqual(b, d, &flg); CHKERRQ(ierr);
1505d4bb536fSBarry Smith   }
1506ca161407SBarry Smith   ierr = MPI_Allreduce(&flg, flag, 1, MPI_INT, MPI_LAND, A->comm);CHKERRQ(ierr);
15073a40ed3dSBarry Smith   PetscFunctionReturn(0);
1508d4bb536fSBarry Smith }
1509d4bb536fSBarry Smith 
15108f6be9afSLois Curfman McInnes extern int MatConvertSameType_MPIAIJ(Mat,Mat *,int);
15112f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int);
15120a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring);
15130a198c4cSBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **);
15146a6a5d1dSBarry Smith extern int MatGetSubMatrix_MPIAIJ (Mat ,IS,IS,int,MatGetSubMatrixCall,Mat *);
151500e6dbe6SBarry Smith 
15168a729477SBarry Smith /* -------------------------------------------------------------------*/
15172ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ,
151839e00950SLois Curfman McInnes        MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ,
151944a69424SLois Curfman McInnes        MatMult_MPIAIJ,MatMultAdd_MPIAIJ,
152044a69424SLois Curfman McInnes        MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ,
152136ce4990SBarry Smith        0,0,
152236ce4990SBarry Smith        0,0,
152336ce4990SBarry Smith        0,0,
152444a69424SLois Curfman McInnes        MatRelax_MPIAIJ,
1525b7c46309SBarry Smith        MatTranspose_MPIAIJ,
1526d4bb536fSBarry Smith        MatGetInfo_MPIAIJ,MatEqual_MPIAIJ,
1527a008b906SSatish Balay        MatGetDiagonal_MPIAIJ,MatDiagonalScale_MPIAIJ,MatNorm_MPIAIJ,
1528ee50ffe9SBarry Smith        MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ,
15291eb62cbbSBarry Smith        0,
1530299609e3SLois Curfman McInnes        MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ,
153136ce4990SBarry Smith        0,0,0,0,
1532d1710a03SLois Curfman McInnes        MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ,
153336ce4990SBarry Smith        0,0,
153494a9d846SBarry Smith        0,0,MatConvertSameType_MPIAIJ,0,0,
1535b49de8d1SLois Curfman McInnes        0,0,0,
1536598137ffSSatish Balay        MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0,
1537052efed2SBarry Smith        MatPrintHelp_MPIAIJ,
15383b2fbd54SBarry Smith        MatScale_MPIAIJ,0,0,0,
15390a198c4cSBarry Smith        MatGetBlockSize_MPIAIJ,0,0,0,0,
154000e6dbe6SBarry Smith        MatFDColoringCreate_MPIAIJ,0,MatSetUnfactored_MPIAIJ,
1541ca161407SBarry Smith        0,0,MatGetSubMatrix_MPIAIJ};
154236ce4990SBarry Smith 
15438a729477SBarry Smith 
15445615d1e5SSatish Balay #undef __FUNC__
15455615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ"
15461987afe7SBarry Smith /*@C
1547ff756334SLois Curfman McInnes    MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format
15483a511b96SLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
15493a511b96SLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameters
15503a511b96SLois Curfman McInnes    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
15513a511b96SLois Curfman McInnes    performance can be increased by more than a factor of 50.
15528a729477SBarry Smith 
1553db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1554db81eaa0SLois Curfman McInnes 
15558a729477SBarry Smith    Input Parameters:
1556db81eaa0SLois Curfman McInnes +  comm - MPI communicator
15577d3e4905SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
155892e8d321SLois Curfman McInnes            This value should be the same as the local size used in creating the
155992e8d321SLois Curfman McInnes            y vector for the matrix-vector product y = Ax.
15601a3896d6SBarry Smith .  n - This value should be the same as the local size used in creating the
15611a3896d6SBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
15621a3896d6SBarry Smith        calculated if N is given) For square matrices n is almost always m.
15637d3e4905SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
156492e8d321SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1565ab693e5aSLois Curfman McInnes .  d_nz - number of nonzeros per row in diagonal portion of local submatrix
1566ff756334SLois Curfman McInnes            (same for all local rows)
15672bd5e0b2SLois Curfman McInnes .  d_nzz - array containing the number of nonzeros in the various rows of the
156892e8d321SLois Curfman McInnes            diagonal portion of the local submatrix (possibly different for each row)
15692bd5e0b2SLois Curfman McInnes            or PETSC_NULL. You must leave room for the diagonal entry even if
15702bd5e0b2SLois Curfman McInnes            it is zero.
15712bd5e0b2SLois Curfman McInnes .  o_nz - number of nonzeros per row in the off-diagonal portion of local
1572ab693e5aSLois Curfman McInnes            submatrix (same for all local rows).
1573db81eaa0SLois Curfman McInnes -  o_nzz - array containing the number of nonzeros in the various rows of the
15742bd5e0b2SLois Curfman McInnes            off-diagonal portion of the local submatrix (possibly different for
15752bd5e0b2SLois Curfman McInnes            each row) or PETSC_NULL.
15768a729477SBarry Smith 
1577ff756334SLois Curfman McInnes    Output Parameter:
157844cd7ae7SLois Curfman McInnes .  A - the matrix
15798a729477SBarry Smith 
1580b259b22eSLois Curfman McInnes    Notes:
1581ff756334SLois Curfman McInnes    The AIJ format (also called the Yale sparse matrix format or
1582ff756334SLois Curfman McInnes    compressed row storage), is fully compatible with standard Fortran 77
15830002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
15840002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
1585ff756334SLois Curfman McInnes 
1586ff756334SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
1587ff756334SLois Curfman McInnes    (possibly both).
1588ff756334SLois Curfman McInnes 
15895511cfe3SLois Curfman McInnes    By default, this format uses inodes (identical nodes) when possible.
15905511cfe3SLois Curfman McInnes    We search for consecutive rows with the same nonzero structure, thereby
15915511cfe3SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
15925511cfe3SLois Curfman McInnes 
15935511cfe3SLois Curfman McInnes    Options Database Keys:
1594db81eaa0SLois Curfman McInnes +  -mat_aij_no_inode  - Do not use inodes
1595db81eaa0SLois Curfman McInnes .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
1596db81eaa0SLois Curfman McInnes -  -mat_aij_oneindex - Internally use indexing starting at 1
1597db81eaa0SLois Curfman McInnes         rather than 0.  Note that when calling MatSetValues(),
1598db81eaa0SLois Curfman McInnes         the user still MUST index entries starting at 0!
15995511cfe3SLois Curfman McInnes 
1600e0245417SLois Curfman McInnes    Storage Information:
1601e0245417SLois Curfman McInnes    For a square global matrix we define each processor's diagonal portion
1602e0245417SLois Curfman McInnes    to be its local rows and the corresponding columns (a square submatrix);
1603e0245417SLois Curfman McInnes    each processor's off-diagonal portion encompasses the remainder of the
1604e0245417SLois Curfman McInnes    local matrix (a rectangular submatrix).
1605e0245417SLois Curfman McInnes 
1606e0245417SLois Curfman McInnes    The user can specify preallocated storage for the diagonal part of
16075ace5be8SLois Curfman McInnes    the local submatrix with either d_nz or d_nnz (not both).  Set
16085ace5be8SLois Curfman McInnes    d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
16095ace5be8SLois Curfman McInnes    memory allocation.  Likewise, specify preallocated storage for the
16105ace5be8SLois Curfman McInnes    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
1611ff756334SLois Curfman McInnes 
16125511cfe3SLois Curfman McInnes    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
16135511cfe3SLois Curfman McInnes    the figure below we depict these three local rows and all columns (0-11).
16142191d07cSBarry Smith 
1615db81eaa0SLois Curfman McInnes .vb
1616db81eaa0SLois Curfman McInnes              0 1 2 3 4 5 6 7 8 9 10 11
1617db81eaa0SLois Curfman McInnes             -------------------
1618db81eaa0SLois Curfman McInnes      row 3  |  o o o d d d o o o o o o
1619db81eaa0SLois Curfman McInnes      row 4  |  o o o d d d o o o o o o
1620db81eaa0SLois Curfman McInnes      row 5  |  o o o d d d o o o o o o
1621db81eaa0SLois Curfman McInnes             -------------------
1622db81eaa0SLois Curfman McInnes .ve
1623b810aeb4SBarry Smith 
16245511cfe3SLois Curfman McInnes    Thus, any entries in the d locations are stored in the d (diagonal)
16255511cfe3SLois Curfman McInnes    submatrix, and any entries in the o locations are stored in the
16265511cfe3SLois Curfman McInnes    o (off-diagonal) submatrix.  Note that the d and the o submatrices are
16275511cfe3SLois Curfman McInnes    stored simply in the MATSEQAIJ format for compressed row storage.
16285511cfe3SLois Curfman McInnes 
16295511cfe3SLois Curfman McInnes    Now d_nz should indicate the number of nonzeros per row in the d matrix,
16305511cfe3SLois Curfman McInnes    and o_nz should indicate the number of nonzeros per row in the o matrix.
16315511cfe3SLois Curfman McInnes    In general, for PDE problems in which most nonzeros are near the diagonal,
16323d323bbdSBarry Smith    one expects d_nz >> o_nz. For large problems you MUST preallocate memory
163392e8d321SLois Curfman McInnes    or you will get TERRIBLE performance; see the users' manual chapter on
16346da5968aSLois Curfman McInnes    matrices.
16353a511b96SLois Curfman McInnes 
1636dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel
1637ff756334SLois Curfman McInnes 
1638fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
16398a729477SBarry Smith @*/
1640e1311b90SBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A)
16418a729477SBarry Smith {
164244cd7ae7SLois Curfman McInnes   Mat          B;
164344cd7ae7SLois Curfman McInnes   Mat_MPIAIJ   *b;
164436ce4990SBarry Smith   int          ierr, i,sum[2],work[2],size;
1645416022c9SBarry Smith 
16463a40ed3dSBarry Smith   PetscFunctionBegin;
16473914022bSBarry Smith   MPI_Comm_size(comm,&size);
16483914022bSBarry Smith   if (size == 1) {
16493914022bSBarry Smith     if (M == PETSC_DECIDE) M = m;
16503914022bSBarry Smith     if (N == PETSC_DECIDE) N = n;
16513914022bSBarry Smith     ierr = MatCreateSeqAIJ(comm,M,N,d_nz,d_nnz,A); CHKERRQ(ierr);
16523a40ed3dSBarry Smith     PetscFunctionReturn(0);
16533914022bSBarry Smith   }
16543914022bSBarry Smith 
165544cd7ae7SLois Curfman McInnes   *A = 0;
1656f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIAIJ,comm,MatDestroy,MatView);
165744cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
165844cd7ae7SLois Curfman McInnes   B->data       = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b);
165944cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_MPIAIJ));
1660f830108cSBarry Smith   PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps));
1661e1311b90SBarry Smith   B->ops->destroy    = MatDestroy_MPIAIJ;
1662e1311b90SBarry Smith   B->ops->view       = MatView_MPIAIJ;
166344cd7ae7SLois Curfman McInnes   B->factor     = 0;
166444cd7ae7SLois Curfman McInnes   B->assembled  = PETSC_FALSE;
166590f02eecSBarry Smith   B->mapping    = 0;
1666d6dfbf8fSBarry Smith 
166747794344SBarry Smith   B->insertmode = NOT_SET_VALUES;
16689eb4d147SSatish Balay   b->size       = size;
166944cd7ae7SLois Curfman McInnes   MPI_Comm_rank(comm,&b->rank);
16701eb62cbbSBarry Smith 
16713a40ed3dSBarry Smith   if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL)) {
1672a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz");
16733a40ed3dSBarry Smith   }
16741987afe7SBarry Smith 
1675dbd7a890SLois Curfman McInnes   if (M == PETSC_DECIDE || N == PETSC_DECIDE) {
16761eb62cbbSBarry Smith     work[0] = m; work[1] = n;
1677ca161407SBarry Smith     ierr = MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm );CHKERRQ(ierr);
1678dbd7a890SLois Curfman McInnes     if (M == PETSC_DECIDE) M = sum[0];
1679dbd7a890SLois Curfman McInnes     if (N == PETSC_DECIDE) N = sum[1];
16801eb62cbbSBarry Smith   }
168144cd7ae7SLois Curfman McInnes   if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);}
168244cd7ae7SLois Curfman McInnes   if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);}
168344cd7ae7SLois Curfman McInnes   b->m = m; B->m = m;
168444cd7ae7SLois Curfman McInnes   b->n = n; B->n = n;
168544cd7ae7SLois Curfman McInnes   b->N = N; B->N = N;
168644cd7ae7SLois Curfman McInnes   b->M = M; B->M = M;
16871eb62cbbSBarry Smith 
16881eb62cbbSBarry Smith   /* build local table of row and column ownerships */
168944cd7ae7SLois Curfman McInnes   b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners);
1690f09e8eb9SSatish Balay   PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ));
1691603f58a4SSatish Balay   b->cowners = b->rowners + b->size + 2;
1692ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
169344cd7ae7SLois Curfman McInnes   b->rowners[0] = 0;
169444cd7ae7SLois Curfman McInnes   for ( i=2; i<=b->size; i++ ) {
169544cd7ae7SLois Curfman McInnes     b->rowners[i] += b->rowners[i-1];
16968a729477SBarry Smith   }
169744cd7ae7SLois Curfman McInnes   b->rstart = b->rowners[b->rank];
169844cd7ae7SLois Curfman McInnes   b->rend   = b->rowners[b->rank+1];
1699ca161407SBarry Smith   ierr = MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
170044cd7ae7SLois Curfman McInnes   b->cowners[0] = 0;
170144cd7ae7SLois Curfman McInnes   for ( i=2; i<=b->size; i++ ) {
170244cd7ae7SLois Curfman McInnes     b->cowners[i] += b->cowners[i-1];
17038a729477SBarry Smith   }
170444cd7ae7SLois Curfman McInnes   b->cstart = b->cowners[b->rank];
170544cd7ae7SLois Curfman McInnes   b->cend   = b->cowners[b->rank+1];
17068a729477SBarry Smith 
17075ace5be8SLois Curfman McInnes   if (d_nz == PETSC_DEFAULT) d_nz = 5;
1708029af93fSBarry Smith   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr);
170944cd7ae7SLois Curfman McInnes   PLogObjectParent(B,b->A);
17107b8455f0SLois Curfman McInnes   if (o_nz == PETSC_DEFAULT) o_nz = 0;
1711029af93fSBarry Smith   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr);
171244cd7ae7SLois Curfman McInnes   PLogObjectParent(B,b->B);
17138a729477SBarry Smith 
17141eb62cbbSBarry Smith   /* build cache for off array entries formed */
171544cd7ae7SLois Curfman McInnes   ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr);
171690f02eecSBarry Smith   b->donotstash  = 0;
171744cd7ae7SLois Curfman McInnes   b->colmap      = 0;
171844cd7ae7SLois Curfman McInnes   b->garray      = 0;
171944cd7ae7SLois Curfman McInnes   b->roworiented = 1;
17208a729477SBarry Smith 
17211eb62cbbSBarry Smith   /* stuff used for matrix vector multiply */
172244cd7ae7SLois Curfman McInnes   b->lvec      = 0;
172344cd7ae7SLois Curfman McInnes   b->Mvctx     = 0;
17248a729477SBarry Smith 
17257a0afa10SBarry Smith   /* stuff for MatGetRow() */
172644cd7ae7SLois Curfman McInnes   b->rowindices   = 0;
172744cd7ae7SLois Curfman McInnes   b->rowvalues    = 0;
172844cd7ae7SLois Curfman McInnes   b->getrowactive = PETSC_FALSE;
17297a0afa10SBarry Smith 
173044cd7ae7SLois Curfman McInnes   *A = B;
17313a40ed3dSBarry Smith   PetscFunctionReturn(0);
1732d6dfbf8fSBarry Smith }
1733c74985f6SBarry Smith 
17345615d1e5SSatish Balay #undef __FUNC__
17355615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIAIJ"
17368f6be9afSLois Curfman McInnes int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues)
1737d6dfbf8fSBarry Smith {
1738d6dfbf8fSBarry Smith   Mat        mat;
1739416022c9SBarry Smith   Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data;
1740a1b97e82SLois Curfman McInnes   int        ierr, len=0, flg;
1741d6dfbf8fSBarry Smith 
17423a40ed3dSBarry Smith   PetscFunctionBegin;
1743416022c9SBarry Smith   *newmat       = 0;
1744f830108cSBarry Smith   PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIAIJ,matin->comm,MatDestroy,MatView);
1745a5a9c739SBarry Smith   PLogObjectCreate(mat);
17460452661fSBarry Smith   mat->data       = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a);
1747f830108cSBarry Smith   PetscMemcpy(mat->ops,&MatOps,sizeof(struct _MatOps));
1748e1311b90SBarry Smith   mat->ops->destroy    = MatDestroy_MPIAIJ;
1749e1311b90SBarry Smith   mat->ops->view       = MatView_MPIAIJ;
1750d6dfbf8fSBarry Smith   mat->factor     = matin->factor;
1751c456f294SBarry Smith   mat->assembled  = PETSC_TRUE;
1752d6dfbf8fSBarry Smith 
175344cd7ae7SLois Curfman McInnes   a->m = mat->m   = oldmat->m;
175444cd7ae7SLois Curfman McInnes   a->n = mat->n   = oldmat->n;
175544cd7ae7SLois Curfman McInnes   a->M = mat->M   = oldmat->M;
175644cd7ae7SLois Curfman McInnes   a->N = mat->N   = oldmat->N;
1757d6dfbf8fSBarry Smith 
1758416022c9SBarry Smith   a->rstart       = oldmat->rstart;
1759416022c9SBarry Smith   a->rend         = oldmat->rend;
1760416022c9SBarry Smith   a->cstart       = oldmat->cstart;
1761416022c9SBarry Smith   a->cend         = oldmat->cend;
176217699dbbSLois Curfman McInnes   a->size         = oldmat->size;
176317699dbbSLois Curfman McInnes   a->rank         = oldmat->rank;
176447794344SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1765bcd2baecSBarry Smith   a->rowvalues    = 0;
1766bcd2baecSBarry Smith   a->getrowactive = PETSC_FALSE;
1767d6dfbf8fSBarry Smith 
1768603f58a4SSatish Balay   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners);
1769f09e8eb9SSatish Balay   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAIJ));
1770603f58a4SSatish Balay   a->cowners = a->rowners + a->size + 2;
1771603f58a4SSatish Balay   PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int));
1772416022c9SBarry Smith   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
17732ee70a88SLois Curfman McInnes   if (oldmat->colmap) {
17740452661fSBarry Smith     a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap);
1775416022c9SBarry Smith     PLogObjectMemory(mat,(a->N)*sizeof(int));
1776416022c9SBarry Smith     PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int));
1777416022c9SBarry Smith   } else a->colmap = 0;
17783f41c07dSBarry Smith   if (oldmat->garray) {
17793f41c07dSBarry Smith     len = ((Mat_SeqAIJ *) (oldmat->B->data))->n;
17803f41c07dSBarry Smith     a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray);
1781464493b3SBarry Smith     PLogObjectMemory(mat,len*sizeof(int));
17823f41c07dSBarry Smith     if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int));
1783416022c9SBarry Smith   } else a->garray = 0;
1784d6dfbf8fSBarry Smith 
1785416022c9SBarry Smith   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
1786416022c9SBarry Smith   PLogObjectParent(mat,a->lvec);
1787a56f8943SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
1788416022c9SBarry Smith   PLogObjectParent(mat,a->Mvctx);
1789416022c9SBarry Smith   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
1790416022c9SBarry Smith   PLogObjectParent(mat,a->A);
1791416022c9SBarry Smith   ierr =  MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr);
1792416022c9SBarry Smith   PLogObjectParent(mat,a->B);
17935dd7a6c7SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
179425cdf11fSBarry Smith   if (flg) {
1795682d7d0cSBarry Smith     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
1796682d7d0cSBarry Smith   }
17978a729477SBarry Smith   *newmat = mat;
17983a40ed3dSBarry Smith   PetscFunctionReturn(0);
17998a729477SBarry Smith }
1800416022c9SBarry Smith 
180177c4ece6SBarry Smith #include "sys.h"
1802416022c9SBarry Smith 
18035615d1e5SSatish Balay #undef __FUNC__
18045615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ"
180519bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat)
1806416022c9SBarry Smith {
1807d65a2f8fSBarry Smith   Mat          A;
1808d65a2f8fSBarry Smith   Scalar       *vals,*svals;
180919bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1810416022c9SBarry Smith   MPI_Status   status;
18113a40ed3dSBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
181217699dbbSLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
1813d65a2f8fSBarry Smith   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
181419bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
1815416022c9SBarry Smith 
18163a40ed3dSBarry Smith   PetscFunctionBegin;
181717699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank);
181817699dbbSLois Curfman McInnes   if (!rank) {
181990ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
18200752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr);
1821a8c6a408SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object");
1822d64ed03dSBarry Smith     if (header[3] < 0) {
1823a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix in special format on disk, cannot load as MPIAIJ");
1824d64ed03dSBarry Smith     }
18256c5fab8fSBarry Smith   }
18266c5fab8fSBarry Smith 
1827d64ed03dSBarry Smith 
1828ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
1829416022c9SBarry Smith   M = header[1]; N = header[2];
1830416022c9SBarry Smith   /* determine ownership of all rows */
183117699dbbSLois Curfman McInnes   m = M/size + ((M % size) > rank);
18320452661fSBarry Smith   rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
1833ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
1834416022c9SBarry Smith   rowners[0] = 0;
183517699dbbSLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
1836416022c9SBarry Smith     rowners[i] += rowners[i-1];
1837416022c9SBarry Smith   }
183817699dbbSLois Curfman McInnes   rstart = rowners[rank];
183917699dbbSLois Curfman McInnes   rend   = rowners[rank+1];
1840416022c9SBarry Smith 
1841416022c9SBarry Smith   /* distribute row lengths to all processors */
18420452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
1843416022c9SBarry Smith   offlens = ourlens + (rend-rstart);
184417699dbbSLois Curfman McInnes   if (!rank) {
18450452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
18460752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
18470452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts);
184817699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1849ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
18500452661fSBarry Smith     PetscFree(sndcounts);
18513a40ed3dSBarry Smith   } else {
1852ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr);
1853416022c9SBarry Smith   }
1854416022c9SBarry Smith 
185517699dbbSLois Curfman McInnes   if (!rank) {
1856416022c9SBarry Smith     /* calculate the number of nonzeros on each processor */
18570452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz);
1858cddf8d76SBarry Smith     PetscMemzero(procsnz,size*sizeof(int));
185917699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) {
1860416022c9SBarry Smith       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
1861416022c9SBarry Smith         procsnz[i] += rowlengths[j];
1862416022c9SBarry Smith       }
1863416022c9SBarry Smith     }
18640452661fSBarry Smith     PetscFree(rowlengths);
1865416022c9SBarry Smith 
1866416022c9SBarry Smith     /* determine max buffer needed and allocate it */
1867416022c9SBarry Smith     maxnz = 0;
186817699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) {
18690452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
1870416022c9SBarry Smith     }
18710452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols);
1872416022c9SBarry Smith 
1873416022c9SBarry Smith     /* read in my part of the matrix column indices  */
1874416022c9SBarry Smith     nz     = procsnz[0];
18750452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
18760752156aSBarry Smith     ierr   = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr);
1877d65a2f8fSBarry Smith 
1878d65a2f8fSBarry Smith     /* read in every one elses and ship off */
187917699dbbSLois Curfman McInnes     for ( i=1; i<size; i++ ) {
1880d65a2f8fSBarry Smith       nz   = procsnz[i];
18810752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr);
1882ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
1883d65a2f8fSBarry Smith     }
18840452661fSBarry Smith     PetscFree(cols);
18853a40ed3dSBarry Smith   } else {
1886416022c9SBarry Smith     /* determine buffer space needed for message */
1887416022c9SBarry Smith     nz = 0;
1888416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
1889416022c9SBarry Smith       nz += ourlens[i];
1890416022c9SBarry Smith     }
18910452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
1892416022c9SBarry Smith 
1893416022c9SBarry Smith     /* receive message of column indices*/
1894ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1895ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
1896a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
1897416022c9SBarry Smith   }
1898416022c9SBarry Smith 
1899416022c9SBarry Smith   /* loop over local rows, determining number of off diagonal entries */
1900cddf8d76SBarry Smith   PetscMemzero(offlens,m*sizeof(int));
1901416022c9SBarry Smith   jj = 0;
1902416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1903416022c9SBarry Smith     for ( j=0; j<ourlens[i]; j++ ) {
1904d65a2f8fSBarry Smith       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
1905416022c9SBarry Smith       jj++;
1906416022c9SBarry Smith     }
1907416022c9SBarry Smith   }
1908d65a2f8fSBarry Smith 
1909d65a2f8fSBarry Smith   /* create our matrix */
1910416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1911416022c9SBarry Smith     ourlens[i] -= offlens[i];
1912416022c9SBarry Smith   }
1913d65a2f8fSBarry Smith   ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr);
1914d65a2f8fSBarry Smith   A = *newmat;
19156d4a8577SBarry Smith   MatSetOption(A,MAT_COLUMNS_SORTED);
1916d65a2f8fSBarry Smith   for ( i=0; i<m; i++ ) {
1917d65a2f8fSBarry Smith     ourlens[i] += offlens[i];
1918d65a2f8fSBarry Smith   }
1919416022c9SBarry Smith 
192017699dbbSLois Curfman McInnes   if (!rank) {
19210452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
1922416022c9SBarry Smith 
1923416022c9SBarry Smith     /* read in my part of the matrix numerical values  */
1924416022c9SBarry Smith     nz = procsnz[0];
19250752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
1926d65a2f8fSBarry Smith 
1927d65a2f8fSBarry Smith     /* insert into matrix */
1928d65a2f8fSBarry Smith     jj      = rstart;
1929d65a2f8fSBarry Smith     smycols = mycols;
1930d65a2f8fSBarry Smith     svals   = vals;
1931d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1932d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1933d65a2f8fSBarry Smith       smycols += ourlens[i];
1934d65a2f8fSBarry Smith       svals   += ourlens[i];
1935d65a2f8fSBarry Smith       jj++;
1936416022c9SBarry Smith     }
1937416022c9SBarry Smith 
1938d65a2f8fSBarry Smith     /* read in other processors and ship out */
193917699dbbSLois Curfman McInnes     for ( i=1; i<size; i++ ) {
1940416022c9SBarry Smith       nz   = procsnz[i];
19410752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr);
1942ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
1943416022c9SBarry Smith     }
19440452661fSBarry Smith     PetscFree(procsnz);
19453a40ed3dSBarry Smith   } else {
1946d65a2f8fSBarry Smith     /* receive numeric values */
19470452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
1948416022c9SBarry Smith 
1949d65a2f8fSBarry Smith     /* receive message of values*/
1950ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1951ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
1952a8c6a408SBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file");
1953d65a2f8fSBarry Smith 
1954d65a2f8fSBarry Smith     /* insert into matrix */
1955d65a2f8fSBarry Smith     jj      = rstart;
1956d65a2f8fSBarry Smith     smycols = mycols;
1957d65a2f8fSBarry Smith     svals   = vals;
1958d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1959d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1960d65a2f8fSBarry Smith       smycols += ourlens[i];
1961d65a2f8fSBarry Smith       svals   += ourlens[i];
1962d65a2f8fSBarry Smith       jj++;
1963d65a2f8fSBarry Smith     }
1964d65a2f8fSBarry Smith   }
19650452661fSBarry Smith   PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners);
1966d65a2f8fSBarry Smith 
19676d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
19686d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
19693a40ed3dSBarry Smith   PetscFunctionReturn(0);
1970416022c9SBarry Smith }
1971a0ff6018SBarry Smith 
197229da9460SBarry Smith #undef __FUNC__
197329da9460SBarry Smith #define __FUNC__ "MatGetSubMatrix_MPIAIJ"
1974a0ff6018SBarry Smith /*
197529da9460SBarry Smith     Not great since it makes two copies of the submatrix, first an SeqAIJ
197629da9460SBarry Smith   in local and then by concatenating the local matrices the end result.
197729da9460SBarry Smith   Writing it directly would be much like MatGetSubMatrices_MPIAIJ()
1978a0ff6018SBarry Smith */
19796a6a5d1dSBarry Smith int MatGetSubMatrix_MPIAIJ(Mat mat,IS isrow,IS iscol,int csize,MatGetSubMatrixCall call,Mat *newmat)
1980a0ff6018SBarry Smith {
198100e6dbe6SBarry Smith   int        ierr, i, m,n,rstart,row,rend,nz,*cwork,size,rank,j;
1982fee21e36SBarry Smith   Mat        *local,M, Mreuse;
198300e6dbe6SBarry Smith   Scalar     *vwork,*aa;
198400e6dbe6SBarry Smith   MPI_Comm   comm = mat->comm;
198500e6dbe6SBarry Smith   Mat_SeqAIJ *aij;
198600e6dbe6SBarry Smith   int        *ii, *jj,nlocal,*dlens,*olens,dlen,olen,jend;
1987a0ff6018SBarry Smith 
1988a0ff6018SBarry Smith   PetscFunctionBegin;
198900e6dbe6SBarry Smith   MPI_Comm_rank(comm,&rank);
199000e6dbe6SBarry Smith   MPI_Comm_size(comm,&size);
199100e6dbe6SBarry Smith 
1992fee21e36SBarry Smith   if (call ==  MAT_REUSE_MATRIX) {
1993fee21e36SBarry Smith     ierr = PetscObjectQuery((PetscObject)*newmat,"SubMatrix",(PetscObject *)&Mreuse);CHKERRQ(ierr);
1994fee21e36SBarry Smith     if (!Mreuse) SETERRQ(1,1,"Submatrix passed in was not used before, cannot reuse");
1995fee21e36SBarry Smith     local = &Mreuse;
1996fee21e36SBarry Smith     ierr  = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_REUSE_MATRIX,&local);CHKERRQ(ierr);
1997fee21e36SBarry Smith   } else {
1998a0ff6018SBarry Smith     ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&local);CHKERRQ(ierr);
1999fee21e36SBarry Smith     Mreuse = *local;
2000fee21e36SBarry Smith     PetscFree(local);
2001fee21e36SBarry Smith   }
2002a0ff6018SBarry Smith 
2003a0ff6018SBarry Smith   /*
2004a0ff6018SBarry Smith       m - number of local rows
2005a0ff6018SBarry Smith       n - number of columns (same on all processors)
2006a0ff6018SBarry Smith       rstart - first row in new global matrix generated
2007a0ff6018SBarry Smith   */
2008fee21e36SBarry Smith   ierr = MatGetSize(Mreuse,&m,&n);CHKERRQ(ierr);
2009a0ff6018SBarry Smith   if (call == MAT_INITIAL_MATRIX) {
2010fee21e36SBarry Smith     aij = (Mat_SeqAIJ *) (Mreuse)->data;
2011a8c6a408SBarry Smith     if (aij->indexshift) SETERRQ(PETSC_ERR_SUP,1,"No support for index shifted matrix");
201200e6dbe6SBarry Smith     ii  = aij->i;
201300e6dbe6SBarry Smith     jj  = aij->j;
201400e6dbe6SBarry Smith 
2015a0ff6018SBarry Smith     /*
201600e6dbe6SBarry Smith         Determine the number of non-zeros in the diagonal and off-diagonal
201700e6dbe6SBarry Smith         portions of the matrix in order to do correct preallocation
2018a0ff6018SBarry Smith     */
201900e6dbe6SBarry Smith 
202000e6dbe6SBarry Smith     /* first get start and end of "diagonal" columns */
20216a6a5d1dSBarry Smith     if (csize == PETSC_DECIDE) {
202200e6dbe6SBarry Smith       nlocal = n/size + ((n % size) > rank);
20236a6a5d1dSBarry Smith     } else {
20246a6a5d1dSBarry Smith       nlocal = csize;
20256a6a5d1dSBarry Smith     }
2026ca161407SBarry Smith     ierr   = MPI_Scan(&nlocal,&rend,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
202700e6dbe6SBarry Smith     rstart = rend - nlocal;
20286a6a5d1dSBarry Smith     if (rank == size - 1 && rend != n) {
20296a6a5d1dSBarry Smith       SETERRQ(1,1,"Local column sizes do not add up to total number of columns");
20306a6a5d1dSBarry Smith     }
203100e6dbe6SBarry Smith 
203200e6dbe6SBarry Smith     /* next, compute all the lengths */
203300e6dbe6SBarry Smith     dlens = (int *) PetscMalloc( (2*m+1)*sizeof(int) );CHKPTRQ(dlens);
203400e6dbe6SBarry Smith     olens = dlens + m;
203500e6dbe6SBarry Smith     for ( i=0; i<m; i++ ) {
203600e6dbe6SBarry Smith       jend = ii[i+1] - ii[i];
203700e6dbe6SBarry Smith       olen = 0;
203800e6dbe6SBarry Smith       dlen = 0;
203900e6dbe6SBarry Smith       for ( j=0; j<jend; j++ ) {
204000e6dbe6SBarry Smith         if ( *jj < rstart || *jj >= rend) olen++;
204100e6dbe6SBarry Smith         else dlen++;
204200e6dbe6SBarry Smith         jj++;
204300e6dbe6SBarry Smith       }
204400e6dbe6SBarry Smith       olens[i] = olen;
204500e6dbe6SBarry Smith       dlens[i] = dlen;
204600e6dbe6SBarry Smith     }
204700e6dbe6SBarry Smith     ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,PETSC_DECIDE,n,0,dlens,0,olens,&M);CHKERRQ(ierr);
204800e6dbe6SBarry Smith     PetscFree(dlens);
2049a0ff6018SBarry Smith   } else {
2050a0ff6018SBarry Smith     int ml,nl;
2051a0ff6018SBarry Smith 
2052a0ff6018SBarry Smith     M = *newmat;
2053a0ff6018SBarry Smith     ierr = MatGetLocalSize(M,&ml,&nl);CHKERRQ(ierr);
2054a8c6a408SBarry Smith     if (ml != m) SETERRQ(PETSC_ERR_ARG_SIZ,1,"Previous matrix must be same size/layout as request");
2055a0ff6018SBarry Smith     ierr = MatZeroEntries(M);CHKERRQ(ierr);
2056c48de900SBarry Smith     /*
2057c48de900SBarry Smith          The next two lines are needed so we may call MatSetValues_MPIAIJ() below directly,
2058c48de900SBarry Smith        rather than the slower MatSetValues().
2059c48de900SBarry Smith     */
2060c48de900SBarry Smith     M->was_assembled = PETSC_TRUE;
2061c48de900SBarry Smith     M->assembled     = PETSC_FALSE;
2062a0ff6018SBarry Smith   }
2063a0ff6018SBarry Smith   ierr = MatGetOwnershipRange(M,&rstart,&rend); CHKERRQ(ierr);
2064fee21e36SBarry Smith   aij = (Mat_SeqAIJ *) (Mreuse)->data;
2065a8c6a408SBarry Smith   if (aij->indexshift) SETERRQ(PETSC_ERR_SUP,1,"No support for index shifted matrix");
206600e6dbe6SBarry Smith   ii  = aij->i;
206700e6dbe6SBarry Smith   jj  = aij->j;
206800e6dbe6SBarry Smith   aa  = aij->a;
2069a0ff6018SBarry Smith   for (i=0; i<m; i++) {
2070a0ff6018SBarry Smith     row   = rstart + i;
207100e6dbe6SBarry Smith     nz    = ii[i+1] - ii[i];
207200e6dbe6SBarry Smith     cwork = jj;     jj += nz;
207300e6dbe6SBarry Smith     vwork = aa;     aa += nz;
20748c638d02SBarry Smith     ierr = MatSetValues_MPIAIJ(M,1,&row,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
2075a0ff6018SBarry Smith   }
2076a0ff6018SBarry Smith 
2077a0ff6018SBarry Smith   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2078a0ff6018SBarry Smith   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2079a0ff6018SBarry Smith   *newmat = M;
2080fee21e36SBarry Smith 
2081fee21e36SBarry Smith   /* save submatrix used in processor for next request */
2082fee21e36SBarry Smith   if (call ==  MAT_INITIAL_MATRIX) {
2083fee21e36SBarry Smith     ierr = PetscObjectCompose((PetscObject)M,"SubMatrix",(PetscObject)Mreuse);CHKERRQ(ierr);
2084fee21e36SBarry Smith     ierr = PetscObjectDereference((PetscObject)Mreuse);CHKERRQ(ierr);
2085fee21e36SBarry Smith   }
2086fee21e36SBarry Smith 
2087a0ff6018SBarry Smith   PetscFunctionReturn(0);
2088a0ff6018SBarry Smith }
2089