xref: /petsc/src/mat/impls/aij/mpi/mpiaij.c (revision 11ebbc714b45ef4372941d536196407bce829dd6)
1cb512458SBarry Smith #ifndef lint
2*11ebbc71SLois Curfman McInnes static char vcid[] = "$Id: mpiaij.c,v 1.203 1997/05/03 22:43:34 curfman Exp curfman $";
3cb512458SBarry Smith #endif
48a729477SBarry Smith 
53369ce9aSBarry Smith #include "pinclude/pviewer.h"
670f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h"
7f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
8d9942c19SSatish Balay #include "src/inline/spops.h"
98a729477SBarry Smith 
109e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column
119e25ed09SBarry Smith number to the local number in the off-diagonal part of the local
129e25ed09SBarry Smith storage of the matrix.  This is done in a non scable way since the
139e25ed09SBarry Smith length of colmap equals the global matrix length.
149e25ed09SBarry Smith */
155615d1e5SSatish Balay #undef __FUNC__
165eea60f9SBarry Smith #define __FUNC__ "CreateColmap_MPIAIJ_Private" /* ADIC Ignore */
170a198c4cSBarry Smith int CreateColmap_MPIAIJ_Private(Mat mat)
189e25ed09SBarry Smith {
1944a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
20ec8511deSBarry Smith   Mat_SeqAIJ *B = (Mat_SeqAIJ*) aij->B->data;
21905e6a2fSBarry Smith   int        n = B->n,i;
22dbb450caSBarry Smith 
23758f045eSSatish Balay   aij->colmap = (int *) PetscMalloc((aij->N+1)*sizeof(int));CHKPTRQ(aij->colmap);
24464493b3SBarry Smith   PLogObjectMemory(mat,aij->N*sizeof(int));
25cddf8d76SBarry Smith   PetscMemzero(aij->colmap,aij->N*sizeof(int));
26905e6a2fSBarry Smith   for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1;
279e25ed09SBarry Smith   return 0;
289e25ed09SBarry Smith }
299e25ed09SBarry Smith 
302493cbb0SBarry Smith extern int DisAssemble_MPIAIJ(Mat);
312493cbb0SBarry Smith 
325615d1e5SSatish Balay #undef __FUNC__
335615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_MPIAIJ"
348f6be9afSLois Curfman McInnes int MatGetRowIJ_MPIAIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int **ja,
353b2fbd54SBarry Smith                            PetscTruth *done)
36299609e3SLois Curfman McInnes {
37299609e3SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
38299609e3SLois Curfman McInnes   int        ierr;
3917699dbbSLois Curfman McInnes   if (aij->size == 1) {
403b2fbd54SBarry Smith     ierr = MatGetRowIJ(aij->A,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
41e3372554SBarry Smith   } else SETERRQ(1,0,"not supported in parallel");
423b2fbd54SBarry Smith   return 0;
433b2fbd54SBarry Smith }
443b2fbd54SBarry Smith 
455615d1e5SSatish Balay #undef __FUNC__
465eea60f9SBarry Smith #define __FUNC__ "MatRestoreRowIJ_MPIAIJ" /* ADIC Ignore */
478f6be9afSLois Curfman McInnes int MatRestoreRowIJ_MPIAIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int **ja,
483b2fbd54SBarry Smith                                PetscTruth *done)
493b2fbd54SBarry Smith {
503b2fbd54SBarry Smith   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
513b2fbd54SBarry Smith   int        ierr;
523b2fbd54SBarry Smith   if (aij->size == 1) {
533b2fbd54SBarry Smith     ierr = MatRestoreRowIJ(aij->A,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
54e3372554SBarry Smith   } else SETERRQ(1,0,"not supported in parallel");
55299609e3SLois Curfman McInnes   return 0;
56299609e3SLois Curfman McInnes }
57299609e3SLois Curfman McInnes 
580520107fSSatish Balay #define CHUNKSIZE   15
5930770e4dSSatish Balay #define MatSetValues_SeqAIJ_A_Private(row,col,value,addv) \
600520107fSSatish Balay { \
610520107fSSatish Balay  \
620520107fSSatish Balay     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift; \
6330770e4dSSatish Balay     rmax = aimax[row]; nrow = ailen[row];  \
64f5e9677aSSatish Balay     col1 = col - shift; \
65f5e9677aSSatish Balay      \
660520107fSSatish Balay       for ( _i=0; _i<nrow; _i++ ) { \
67f5e9677aSSatish Balay         if (rp[_i] > col1) break; \
68f5e9677aSSatish Balay         if (rp[_i] == col1) { \
690520107fSSatish Balay           if (addv == ADD_VALUES) ap[_i] += value;   \
700520107fSSatish Balay           else                  ap[_i] = value; \
7130770e4dSSatish Balay           goto a_noinsert; \
720520107fSSatish Balay         } \
730520107fSSatish Balay       }  \
7489280ab3SLois Curfman McInnes       if (nonew == 1) goto a_noinsert; \
75*11ebbc71SLois Curfman McInnes       else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \
760520107fSSatish Balay       if (nrow >= rmax) { \
770520107fSSatish Balay         /* there is no extra room in row, therefore enlarge */ \
780520107fSSatish Balay         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; \
790520107fSSatish Balay         Scalar *new_a; \
800520107fSSatish Balay  \
81*11ebbc71SLois Curfman McInnes         if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \
8289280ab3SLois Curfman McInnes  \
830520107fSSatish Balay         /* malloc new storage space */ \
840520107fSSatish Balay         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); \
850520107fSSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \
860520107fSSatish Balay         new_j   = (int *) (new_a + new_nz); \
870520107fSSatish Balay         new_i   = new_j + new_nz; \
880520107fSSatish Balay  \
890520107fSSatish Balay         /* copy over old data into new slots */ \
900520107fSSatish Balay         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} \
910520107fSSatish Balay         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} \
920520107fSSatish Balay         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); \
930520107fSSatish Balay         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); \
940520107fSSatish Balay         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, \
950520107fSSatish Balay                                                            len*sizeof(int)); \
960520107fSSatish Balay         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); \
970520107fSSatish Balay         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, \
980520107fSSatish Balay                                                            len*sizeof(Scalar));  \
990520107fSSatish Balay         /* free up old matrix storage */ \
100f5e9677aSSatish Balay  \
1010520107fSSatish Balay         PetscFree(a->a);  \
1020520107fSSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} \
1030520107fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;  \
1040520107fSSatish Balay         a->singlemalloc = 1; \
1050520107fSSatish Balay  \
1060520107fSSatish Balay         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift; \
10730770e4dSSatish Balay         rmax = aimax[row] = aimax[row] + CHUNKSIZE; \
1080520107fSSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \
1090520107fSSatish Balay         a->maxnz += CHUNKSIZE; \
1100520107fSSatish Balay         a->reallocs++; \
1110520107fSSatish Balay       } \
1120520107fSSatish Balay       N = nrow++ - 1; a->nz++; \
1130520107fSSatish Balay       /* shift up all the later entries in this row */ \
1140520107fSSatish Balay       for ( ii=N; ii>=_i; ii-- ) { \
1150520107fSSatish Balay         rp[ii+1] = rp[ii]; \
1160520107fSSatish Balay         ap[ii+1] = ap[ii]; \
1170520107fSSatish Balay       } \
118f5e9677aSSatish Balay       rp[_i] = col1;  \
1190520107fSSatish Balay       ap[_i] = value;  \
12030770e4dSSatish Balay       a_noinsert: ; \
1210520107fSSatish Balay       ailen[row] = nrow; \
1220520107fSSatish Balay }
1230a198c4cSBarry Smith 
12430770e4dSSatish Balay #define MatSetValues_SeqAIJ_B_Private(row,col,value,addv) \
12530770e4dSSatish Balay { \
12630770e4dSSatish Balay  \
12730770e4dSSatish Balay     rp   = bj + bi[row] + shift; ap = ba + bi[row] + shift; \
12830770e4dSSatish Balay     rmax = bimax[row]; nrow = bilen[row];  \
12930770e4dSSatish Balay     col1 = col - shift; \
13030770e4dSSatish Balay      \
13130770e4dSSatish Balay       for ( _i=0; _i<nrow; _i++ ) { \
13230770e4dSSatish Balay         if (rp[_i] > col1) break; \
13330770e4dSSatish Balay         if (rp[_i] == col1) { \
13430770e4dSSatish Balay           if (addv == ADD_VALUES) ap[_i] += value;   \
13530770e4dSSatish Balay           else                  ap[_i] = value; \
13630770e4dSSatish Balay           goto b_noinsert; \
13730770e4dSSatish Balay         } \
13830770e4dSSatish Balay       }  \
13989280ab3SLois Curfman McInnes       if (nonew == 1) goto b_noinsert; \
140*11ebbc71SLois Curfman McInnes       else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \
14130770e4dSSatish Balay       if (nrow >= rmax) { \
14230770e4dSSatish Balay         /* there is no extra room in row, therefore enlarge */ \
14374c639caSSatish Balay         int    new_nz = bi[b->m] + CHUNKSIZE,len,*new_i,*new_j; \
14430770e4dSSatish Balay         Scalar *new_a; \
14530770e4dSSatish Balay  \
146*11ebbc71SLois Curfman McInnes         if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); \
14789280ab3SLois Curfman McInnes  \
14830770e4dSSatish Balay         /* malloc new storage space */ \
14974c639caSSatish Balay         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(b->m+1)*sizeof(int); \
15030770e4dSSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); \
15130770e4dSSatish Balay         new_j   = (int *) (new_a + new_nz); \
15230770e4dSSatish Balay         new_i   = new_j + new_nz; \
15330770e4dSSatish Balay  \
15430770e4dSSatish Balay         /* copy over old data into new slots */ \
15530770e4dSSatish Balay         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = bi[ii];} \
15674c639caSSatish Balay         for ( ii=row+1; ii<b->m+1; ii++ ) {new_i[ii] = bi[ii]+CHUNKSIZE;} \
15730770e4dSSatish Balay         PetscMemcpy(new_j,bj,(bi[row]+nrow+shift)*sizeof(int)); \
15830770e4dSSatish Balay         len = (new_nz - CHUNKSIZE - bi[row] - nrow - shift); \
15930770e4dSSatish Balay         PetscMemcpy(new_j+bi[row]+shift+nrow+CHUNKSIZE,bj+bi[row]+shift+nrow, \
16030770e4dSSatish Balay                                                            len*sizeof(int)); \
16130770e4dSSatish Balay         PetscMemcpy(new_a,ba,(bi[row]+nrow+shift)*sizeof(Scalar)); \
16230770e4dSSatish Balay         PetscMemcpy(new_a+bi[row]+shift+nrow+CHUNKSIZE,ba+bi[row]+shift+nrow, \
16330770e4dSSatish Balay                                                            len*sizeof(Scalar));  \
16430770e4dSSatish Balay         /* free up old matrix storage */ \
16530770e4dSSatish Balay  \
16674c639caSSatish Balay         PetscFree(b->a);  \
16774c639caSSatish Balay         if (!b->singlemalloc) {PetscFree(b->i);PetscFree(b->j);} \
16874c639caSSatish Balay         ba = b->a = new_a; bi = b->i = new_i; bj = b->j = new_j;  \
16974c639caSSatish Balay         b->singlemalloc = 1; \
17030770e4dSSatish Balay  \
17130770e4dSSatish Balay         rp   = bj + bi[row] + shift; ap = ba + bi[row] + shift; \
17230770e4dSSatish Balay         rmax = bimax[row] = bimax[row] + CHUNKSIZE; \
17374c639caSSatish Balay         PLogObjectMemory(B,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); \
17474c639caSSatish Balay         b->maxnz += CHUNKSIZE; \
17574c639caSSatish Balay         b->reallocs++; \
17630770e4dSSatish Balay       } \
17774c639caSSatish Balay       N = nrow++ - 1; b->nz++; \
17830770e4dSSatish Balay       /* shift up all the later entries in this row */ \
17930770e4dSSatish Balay       for ( ii=N; ii>=_i; ii-- ) { \
18030770e4dSSatish Balay         rp[ii+1] = rp[ii]; \
18130770e4dSSatish Balay         ap[ii+1] = ap[ii]; \
18230770e4dSSatish Balay       } \
18330770e4dSSatish Balay       rp[_i] = col1;  \
18430770e4dSSatish Balay       ap[_i] = value;  \
18530770e4dSSatish Balay       b_noinsert: ; \
18630770e4dSSatish Balay       bilen[row] = nrow; \
18730770e4dSSatish Balay }
18830770e4dSSatish Balay 
1890520107fSSatish Balay extern int MatSetValues_SeqAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode);
1905615d1e5SSatish Balay #undef __FUNC__
1915615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIAIJ"
1928f6be9afSLois Curfman McInnes int MatSetValues_MPIAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv)
1938a729477SBarry Smith {
19444a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1954b0e389bSBarry Smith   Scalar     value;
1961eb62cbbSBarry Smith   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
1971eb62cbbSBarry Smith   int        cstart = aij->cstart, cend = aij->cend,row,col;
198905e6a2fSBarry Smith   int        roworiented = aij->roworiented;
1998a729477SBarry Smith 
2000520107fSSatish Balay   /* Some Variables required in the macro */
2014ee7247eSSatish Balay   Mat        A = aij->A;
2024ee7247eSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
20330770e4dSSatish Balay   int        *aimax = a->imax, *ai = a->i, *ailen = a->ilen,*aj = a->j;
20430770e4dSSatish Balay   Scalar     *aa = a->a;
20530770e4dSSatish Balay 
20630770e4dSSatish Balay   Mat        B = aij->B;
20730770e4dSSatish Balay   Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data;
20830770e4dSSatish Balay   int        *bimax = b->imax, *bi = b->i, *bilen = b->ilen,*bj = b->j;
20930770e4dSSatish Balay   Scalar     *ba = b->a;
21030770e4dSSatish Balay 
2114ee7247eSSatish Balay   int        *rp,ii,nrow,_i,rmax, N, col1;
21230770e4dSSatish Balay   int        nonew = a->nonew,shift = a->indexshift;
21330770e4dSSatish Balay   Scalar     *ap;
2144ee7247eSSatish Balay 
2158a729477SBarry Smith   for ( i=0; i<m; i++ ) {
2160a198c4cSBarry Smith #if defined(PETSC_BOPT_g)
217e3372554SBarry Smith     if (im[i] < 0) SETERRQ(1,0,"Negative row");
218e3372554SBarry Smith     if (im[i] >= aij->M) SETERRQ(1,0,"Row too large");
2190a198c4cSBarry Smith #endif
2204b0e389bSBarry Smith     if (im[i] >= rstart && im[i] < rend) {
2214b0e389bSBarry Smith       row = im[i] - rstart;
2221eb62cbbSBarry Smith       for ( j=0; j<n; j++ ) {
2234b0e389bSBarry Smith         if (in[j] >= cstart && in[j] < cend){
2244b0e389bSBarry Smith           col = in[j] - cstart;
2254b0e389bSBarry Smith           if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
22630770e4dSSatish Balay           MatSetValues_SeqAIJ_A_Private(row,col,value,addv);
2270520107fSSatish Balay           /* ierr = MatSetValues_SeqAIJ(aij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */
2281eb62cbbSBarry Smith         }
2290a198c4cSBarry Smith #if defined(PETSC_BOPT_g)
230e3372554SBarry Smith         else if (in[j] < 0) {SETERRQ(1,0,"Negative column");}
231e3372554SBarry Smith         else if (in[j] >= aij->N) {SETERRQ(1,0,"Col too large");}
2320a198c4cSBarry Smith #endif
2331eb62cbbSBarry Smith         else {
234227d817aSBarry Smith           if (mat->was_assembled) {
235905e6a2fSBarry Smith             if (!aij->colmap) {
236905e6a2fSBarry Smith               ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
237905e6a2fSBarry Smith             }
238905e6a2fSBarry Smith             col = aij->colmap[in[j]] - 1;
239ec8511deSBarry Smith             if (col < 0 && !((Mat_SeqAIJ*)(aij->A->data))->nonew) {
2402493cbb0SBarry Smith               ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
2414b0e389bSBarry Smith               col =  in[j];
242d6dfbf8fSBarry Smith             }
2439e25ed09SBarry Smith           }
2444b0e389bSBarry Smith           else col = in[j];
2454b0e389bSBarry Smith           if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
24630770e4dSSatish Balay           MatSetValues_SeqAIJ_B_Private(row,col,value,addv);
24730770e4dSSatish Balay           /* ierr = MatSetValues_SeqAIJ(aij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */
2481eb62cbbSBarry Smith         }
2491eb62cbbSBarry Smith       }
2501eb62cbbSBarry Smith     }
2511eb62cbbSBarry Smith     else {
25290f02eecSBarry Smith       if (roworiented && !aij->donotstash) {
2534b0e389bSBarry Smith         ierr = StashValues_Private(&aij->stash,im[i],n,in,v+i*n,addv);CHKERRQ(ierr);
2544b0e389bSBarry Smith       }
2554b0e389bSBarry Smith       else {
25690f02eecSBarry Smith         if (!aij->donotstash) {
2574b0e389bSBarry Smith           row = im[i];
2584b0e389bSBarry Smith           for ( j=0; j<n; j++ ) {
2594b0e389bSBarry Smith             ierr = StashValues_Private(&aij->stash,row,1,in+j,v+i+j*m,addv);CHKERRQ(ierr);
2604b0e389bSBarry Smith           }
2614b0e389bSBarry Smith         }
2621eb62cbbSBarry Smith       }
2638a729477SBarry Smith     }
26490f02eecSBarry Smith   }
2658a729477SBarry Smith   return 0;
2668a729477SBarry Smith }
2678a729477SBarry Smith 
2685615d1e5SSatish Balay #undef __FUNC__
2695615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIAIJ"
2708f6be9afSLois Curfman McInnes int MatGetValues_MPIAIJ(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
271b49de8d1SLois Curfman McInnes {
272b49de8d1SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
273b49de8d1SLois Curfman McInnes   int        ierr,i,j, rstart = aij->rstart, rend = aij->rend;
274b49de8d1SLois Curfman McInnes   int        cstart = aij->cstart, cend = aij->cend,row,col;
275b49de8d1SLois Curfman McInnes 
276b49de8d1SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
277e3372554SBarry Smith     if (idxm[i] < 0) SETERRQ(1,0,"Negative row");
278e3372554SBarry Smith     if (idxm[i] >= aij->M) SETERRQ(1,0,"Row too large");
279b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
280b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
281b49de8d1SLois Curfman McInnes       for ( j=0; j<n; j++ ) {
282e3372554SBarry Smith         if (idxn[j] < 0) SETERRQ(1,0,"Negative column");
283e3372554SBarry Smith         if (idxn[j] >= aij->N) SETERRQ(1,0,"Col too large");
284b49de8d1SLois Curfman McInnes         if (idxn[j] >= cstart && idxn[j] < cend){
285b49de8d1SLois Curfman McInnes           col = idxn[j] - cstart;
286b49de8d1SLois Curfman McInnes           ierr = MatGetValues(aij->A,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr);
287b49de8d1SLois Curfman McInnes         }
288b49de8d1SLois Curfman McInnes         else {
289905e6a2fSBarry Smith           if (!aij->colmap) {
290905e6a2fSBarry Smith             ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
291905e6a2fSBarry Smith           }
292905e6a2fSBarry Smith           col = aij->colmap[idxn[j]] - 1;
293e60e1c95SSatish Balay           if ((col < 0) || (aij->garray[col] != idxn[j])) *(v+i*n+j) = 0.0;
294d9d09a02SSatish Balay           else {
295b49de8d1SLois Curfman McInnes             ierr = MatGetValues(aij->B,1,&row,1,&col,v+i*n+j); CHKERRQ(ierr);
296b49de8d1SLois Curfman McInnes           }
297b49de8d1SLois Curfman McInnes         }
298b49de8d1SLois Curfman McInnes       }
299d9d09a02SSatish Balay     }
300b49de8d1SLois Curfman McInnes     else {
301e3372554SBarry Smith       SETERRQ(1,0,"Only local values currently supported");
302b49de8d1SLois Curfman McInnes     }
303b49de8d1SLois Curfman McInnes   }
304b49de8d1SLois Curfman McInnes   return 0;
305b49de8d1SLois Curfman McInnes }
306b49de8d1SLois Curfman McInnes 
3075615d1e5SSatish Balay #undef __FUNC__
3085615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIAIJ"
3098f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIAIJ(Mat mat,MatAssemblyType mode)
3108a729477SBarry Smith {
31144a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
312d6dfbf8fSBarry Smith   MPI_Comm    comm = mat->comm;
31317699dbbSLois Curfman McInnes   int         size = aij->size, *owners = aij->rowners;
31417699dbbSLois Curfman McInnes   int         rank = aij->rank,tag = mat->tag, *owner,*starts,count,ierr;
3151eb62cbbSBarry Smith   MPI_Request *send_waits,*recv_waits;
3166abc6512SBarry Smith   int         *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work;
3171eb62cbbSBarry Smith   InsertMode  addv;
3181eb62cbbSBarry Smith   Scalar      *rvalues,*svalues;
3191eb62cbbSBarry Smith 
3201eb62cbbSBarry Smith   /* make sure all processors are either in INSERTMODE or ADDMODE */
32147794344SBarry Smith   MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
322dbb450caSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
323e3372554SBarry Smith     SETERRQ(1,0,"Some processors inserted others added");
3241eb62cbbSBarry Smith   }
32547794344SBarry Smith   mat->insertmode = addv; /* in case this processor had no cache */
3261eb62cbbSBarry Smith 
3271eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
3280452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
329cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
3300452661fSBarry Smith   owner = (int *) PetscMalloc( (aij->stash.n+1)*sizeof(int) ); CHKPTRQ(owner);
3311eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
3321eb62cbbSBarry Smith     idx = aij->stash.idx[i];
33317699dbbSLois Curfman McInnes     for ( j=0; j<size; j++ ) {
3341eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
3351eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; break;
3368a729477SBarry Smith       }
3378a729477SBarry Smith     }
3388a729477SBarry Smith   }
33917699dbbSLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
3401eb62cbbSBarry Smith 
3411eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
3420452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
34317699dbbSLois Curfman McInnes   MPI_Allreduce(procs, work,size,MPI_INT,MPI_SUM,comm);
34417699dbbSLois Curfman McInnes   nreceives = work[rank];
34517699dbbSLois Curfman McInnes   MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);
34617699dbbSLois Curfman McInnes   nmax = work[rank];
3470452661fSBarry Smith   PetscFree(work);
3481eb62cbbSBarry Smith 
3491eb62cbbSBarry Smith   /* post receives:
3501eb62cbbSBarry Smith        1) each message will consist of ordered pairs
3511eb62cbbSBarry Smith      (global index,value) we store the global index as a double
352d6dfbf8fSBarry Smith      to simplify the message passing.
3531eb62cbbSBarry Smith        2) since we don't know how long each individual message is we
3541eb62cbbSBarry Smith      allocate the largest needed buffer for each receive. Potentially
3551eb62cbbSBarry Smith      this is a lot of wasted space.
3561eb62cbbSBarry Smith 
3571eb62cbbSBarry Smith 
3581eb62cbbSBarry Smith        This could be done better.
3591eb62cbbSBarry Smith   */
3600452661fSBarry Smith   rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));
36178b31e54SBarry Smith   CHKPTRQ(rvalues);
3620452661fSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));
36378b31e54SBarry Smith   CHKPTRQ(recv_waits);
3641eb62cbbSBarry Smith   for ( i=0; i<nreceives; i++ ) {
365416022c9SBarry Smith     MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,
3661eb62cbbSBarry Smith               comm,recv_waits+i);
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);
3740452661fSBarry Smith   send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));
37578b31e54SBarry Smith   CHKPTRQ(send_waits);
3760452661fSBarry Smith   starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts);
3771eb62cbbSBarry Smith   starts[0] = 0;
37817699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3791eb62cbbSBarry Smith   for ( i=0; i<aij->stash.n; i++ ) {
3801eb62cbbSBarry Smith     svalues[3*starts[owner[i]]]       = (Scalar)  aij->stash.idx[i];
3811eb62cbbSBarry Smith     svalues[3*starts[owner[i]]+1]     = (Scalar)  aij->stash.idy[i];
3821eb62cbbSBarry Smith     svalues[3*(starts[owner[i]]++)+2] =  aij->stash.array[i];
3831eb62cbbSBarry Smith   }
3840452661fSBarry Smith   PetscFree(owner);
3851eb62cbbSBarry Smith   starts[0] = 0;
38617699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
3871eb62cbbSBarry Smith   count = 0;
38817699dbbSLois Curfman McInnes   for ( i=0; i<size; i++ ) {
3891eb62cbbSBarry Smith     if (procs[i]) {
390416022c9SBarry Smith       MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,
3911eb62cbbSBarry Smith                 comm,send_waits+count++);
3921eb62cbbSBarry Smith     }
3931eb62cbbSBarry Smith   }
3940452661fSBarry Smith   PetscFree(starts); PetscFree(nprocs);
3951eb62cbbSBarry Smith 
3961eb62cbbSBarry Smith   /* Free cache space */
39755908cccSLois Curfman McInnes   PLogInfo(mat,"MatAssemblyBegin_MPIAIJ:Number of off-processor values %d\n",aij->stash.n);
39878b31e54SBarry Smith   ierr = StashDestroy_Private(&aij->stash); CHKERRQ(ierr);
3991eb62cbbSBarry Smith 
4001eb62cbbSBarry Smith   aij->svalues    = svalues;    aij->rvalues    = rvalues;
4011eb62cbbSBarry Smith   aij->nsends     = nsends;     aij->nrecvs     = nreceives;
4021eb62cbbSBarry Smith   aij->send_waits = send_waits; aij->recv_waits = recv_waits;
4031eb62cbbSBarry Smith   aij->rmax       = nmax;
4041eb62cbbSBarry Smith 
4051eb62cbbSBarry Smith   return 0;
4061eb62cbbSBarry Smith }
40744a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat);
4081eb62cbbSBarry Smith 
4095615d1e5SSatish Balay #undef __FUNC__
4105615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIAIJ"
4118f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIAIJ(Mat mat,MatAssemblyType mode)
4121eb62cbbSBarry Smith {
41344a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
4141eb62cbbSBarry Smith   MPI_Status  *send_status,recv_status;
415416022c9SBarry Smith   int         imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n, ierr;
416905e6a2fSBarry Smith   int         row,col,other_disassembled;
4171eb62cbbSBarry Smith   Scalar      *values,val;
41847794344SBarry Smith   InsertMode  addv = mat->insertmode;
4191eb62cbbSBarry Smith 
4201eb62cbbSBarry Smith   /*  wait on receives */
4211eb62cbbSBarry Smith   while (count) {
422d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status);
4231eb62cbbSBarry Smith     /* unpack receives into our local space */
424d6dfbf8fSBarry Smith     values = aij->rvalues + 3*imdex*aij->rmax;
425c60448a5SBarry Smith     MPI_Get_count(&recv_status,MPIU_SCALAR,&n);
4261eb62cbbSBarry Smith     n = n/3;
4271eb62cbbSBarry Smith     for ( i=0; i<n; i++ ) {
428227d817aSBarry Smith       row = (int) PetscReal(values[3*i]) - aij->rstart;
429227d817aSBarry Smith       col = (int) PetscReal(values[3*i+1]);
4301eb62cbbSBarry Smith       val = values[3*i+2];
4311eb62cbbSBarry Smith       if (col >= aij->cstart && col < aij->cend) {
4321eb62cbbSBarry Smith         col -= aij->cstart;
4331eb62cbbSBarry Smith         MatSetValues(aij->A,1,&row,1,&col,&val,addv);
4341eb62cbbSBarry Smith       }
4351eb62cbbSBarry 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         }
4461eb62cbbSBarry Smith         MatSetValues(aij->B,1,&row,1,&col,&val,addv);
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);
4561eb62cbbSBarry Smith     MPI_Waitall(aij->nsends,aij->send_waits,send_status);
4570452661fSBarry Smith     PetscFree(send_status);
4581eb62cbbSBarry Smith   }
4590452661fSBarry Smith   PetscFree(aij->send_waits); PetscFree(aij->svalues);
4601eb62cbbSBarry Smith 
46178b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->A,mode); CHKERRQ(ierr);
46278b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->A,mode); CHKERRQ(ierr);
4631eb62cbbSBarry Smith 
4642493cbb0SBarry Smith   /* determine if any processor has disassembled, if so we must
4652493cbb0SBarry Smith      also disassemble ourselfs, in order that we may reassemble. */
466227d817aSBarry Smith   MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);
467227d817aSBarry Smith   if (mat->was_assembled && !other_disassembled) {
4682493cbb0SBarry Smith     ierr = DisAssemble_MPIAIJ(mat); CHKERRQ(ierr);
4692493cbb0SBarry Smith   }
4702493cbb0SBarry Smith 
4716d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
47278b31e54SBarry Smith     ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERRQ(ierr);
4735e42470aSBarry Smith   }
47478b31e54SBarry Smith   ierr = MatAssemblyBegin(aij->B,mode); CHKERRQ(ierr);
47578b31e54SBarry Smith   ierr = MatAssemblyEnd(aij->B,mode); CHKERRQ(ierr);
4765e42470aSBarry Smith 
4777a0afa10SBarry Smith   if (aij->rowvalues) {PetscFree(aij->rowvalues); aij->rowvalues = 0;}
4788a729477SBarry Smith   return 0;
4798a729477SBarry Smith }
4808a729477SBarry Smith 
4815615d1e5SSatish Balay #undef __FUNC__
4825615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIAIJ"
4838f6be9afSLois Curfman McInnes int MatZeroEntries_MPIAIJ(Mat A)
4841eb62cbbSBarry Smith {
48544a69424SLois Curfman McInnes   Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data;
486dbd7a890SLois Curfman McInnes   int        ierr;
48778b31e54SBarry Smith   ierr = MatZeroEntries(l->A); CHKERRQ(ierr);
48878b31e54SBarry Smith   ierr = MatZeroEntries(l->B); CHKERRQ(ierr);
4891eb62cbbSBarry Smith   return 0;
4901eb62cbbSBarry Smith }
4911eb62cbbSBarry Smith 
4921eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the
4931eb62cbbSBarry Smith    matrix is square and the column and row owerships are identical.
4941eb62cbbSBarry Smith    This is a BUG. The only way to fix it seems to be to access
4951eb62cbbSBarry Smith    aij->A and aij->B directly and not through the MatZeroRows()
4961eb62cbbSBarry Smith    routine.
4971eb62cbbSBarry Smith */
4985615d1e5SSatish Balay #undef __FUNC__
4995615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIAIJ"
5008f6be9afSLois Curfman McInnes int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag)
5011eb62cbbSBarry Smith {
50244a69424SLois Curfman McInnes   Mat_MPIAIJ     *l = (Mat_MPIAIJ *) A->data;
50317699dbbSLois Curfman McInnes   int            i,ierr,N, *rows,*owners = l->rowners,size = l->size;
5046abc6512SBarry Smith   int            *procs,*nprocs,j,found,idx,nsends,*work;
50517699dbbSLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
5065392566eSBarry Smith   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
507d6dfbf8fSBarry Smith   int            *lens,imdex,*lrows,*values;
508d6dfbf8fSBarry Smith   MPI_Comm       comm = A->comm;
5091eb62cbbSBarry Smith   MPI_Request    *send_waits,*recv_waits;
5101eb62cbbSBarry Smith   MPI_Status     recv_status,*send_status;
5111eb62cbbSBarry Smith   IS             istmp;
5121eb62cbbSBarry Smith 
51377c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
51478b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
5151eb62cbbSBarry Smith 
5161eb62cbbSBarry Smith   /*  first count number of contributors to each processor */
5170452661fSBarry Smith   nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs);
518cddf8d76SBarry Smith   PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size;
5190452661fSBarry Smith   owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/
5201eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
5211eb62cbbSBarry Smith     idx = rows[i];
5221eb62cbbSBarry Smith     found = 0;
52317699dbbSLois Curfman McInnes     for ( j=0; j<size; j++ ) {
5241eb62cbbSBarry Smith       if (idx >= owners[j] && idx < owners[j+1]) {
5251eb62cbbSBarry Smith         nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break;
5261eb62cbbSBarry Smith       }
5271eb62cbbSBarry Smith     }
528e3372554SBarry Smith     if (!found) SETERRQ(1,0,"Index out of range");
5291eb62cbbSBarry Smith   }
53017699dbbSLois Curfman McInnes   nsends = 0;  for ( i=0; i<size; i++ ) { nsends += procs[i];}
5311eb62cbbSBarry Smith 
5321eb62cbbSBarry Smith   /* inform other processors of number of messages and max length*/
5330452661fSBarry Smith   work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work);
53417699dbbSLois Curfman McInnes   MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);
53517699dbbSLois Curfman McInnes   nrecvs = work[rank];
53617699dbbSLois Curfman McInnes   MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);
53717699dbbSLois Curfman McInnes   nmax = work[rank];
5380452661fSBarry Smith   PetscFree(work);
5391eb62cbbSBarry Smith 
5401eb62cbbSBarry Smith   /* post receives:   */
5410452661fSBarry Smith   rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */
54278b31e54SBarry Smith   CHKPTRQ(rvalues);
5430452661fSBarry Smith   recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));
54478b31e54SBarry Smith   CHKPTRQ(recv_waits);
5451eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
546416022c9SBarry Smith     MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
5471eb62cbbSBarry Smith   }
5481eb62cbbSBarry Smith 
5491eb62cbbSBarry Smith   /* do sends:
5501eb62cbbSBarry Smith       1) starts[i] gives the starting index in svalues for stuff going to
5511eb62cbbSBarry Smith          the ith processor
5521eb62cbbSBarry Smith   */
5530452661fSBarry Smith   svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues);
5540452661fSBarry Smith   send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request));
55578b31e54SBarry Smith   CHKPTRQ(send_waits);
5560452661fSBarry Smith   starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts);
5571eb62cbbSBarry Smith   starts[0] = 0;
55817699dbbSLois Curfman McInnes   for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
5591eb62cbbSBarry Smith   for ( i=0; i<N; i++ ) {
5601eb62cbbSBarry Smith     svalues[starts[owner[i]]++] = rows[i];
5611eb62cbbSBarry Smith   }
5621eb62cbbSBarry Smith   ISRestoreIndices(is,&rows);
5631eb62cbbSBarry Smith 
5641eb62cbbSBarry Smith   starts[0] = 0;
56517699dbbSLois Curfman McInnes   for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];}
5661eb62cbbSBarry Smith   count = 0;
56717699dbbSLois Curfman McInnes   for ( i=0; i<size; i++ ) {
5681eb62cbbSBarry Smith     if (procs[i]) {
569416022c9SBarry Smith       MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);
5701eb62cbbSBarry Smith     }
5711eb62cbbSBarry Smith   }
5720452661fSBarry Smith   PetscFree(starts);
5731eb62cbbSBarry Smith 
57417699dbbSLois Curfman McInnes   base = owners[rank];
5751eb62cbbSBarry Smith 
5761eb62cbbSBarry Smith   /*  wait on receives */
5770452661fSBarry Smith   lens   = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens);
5781eb62cbbSBarry Smith   source = lens + nrecvs;
5791eb62cbbSBarry Smith   count  = nrecvs; slen = 0;
5801eb62cbbSBarry Smith   while (count) {
581d6dfbf8fSBarry Smith     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
5821eb62cbbSBarry Smith     /* unpack receives into our local space */
5831eb62cbbSBarry Smith     MPI_Get_count(&recv_status,MPI_INT,&n);
584d6dfbf8fSBarry Smith     source[imdex]  = recv_status.MPI_SOURCE;
585d6dfbf8fSBarry Smith     lens[imdex]  = n;
5861eb62cbbSBarry Smith     slen += n;
5871eb62cbbSBarry Smith     count--;
5881eb62cbbSBarry Smith   }
5890452661fSBarry Smith   PetscFree(recv_waits);
5901eb62cbbSBarry Smith 
5911eb62cbbSBarry Smith   /* move the data into the send scatter */
5920452661fSBarry Smith   lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows);
5931eb62cbbSBarry Smith   count = 0;
5941eb62cbbSBarry Smith   for ( i=0; i<nrecvs; i++ ) {
5951eb62cbbSBarry Smith     values = rvalues + i*nmax;
5961eb62cbbSBarry Smith     for ( j=0; j<lens[i]; j++ ) {
5971eb62cbbSBarry Smith       lrows[count++] = values[j] - base;
5981eb62cbbSBarry Smith     }
5991eb62cbbSBarry Smith   }
6000452661fSBarry Smith   PetscFree(rvalues); PetscFree(lens);
6010452661fSBarry Smith   PetscFree(owner); PetscFree(nprocs);
6021eb62cbbSBarry Smith 
6031eb62cbbSBarry Smith   /* actually zap the local rows */
604029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
605464493b3SBarry Smith   PLogObjectParent(A,istmp);
6060452661fSBarry Smith   PetscFree(lrows);
60778b31e54SBarry Smith   ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr);
60878b31e54SBarry Smith   ierr = MatZeroRows(l->B,istmp,0); CHKERRQ(ierr);
60978b31e54SBarry Smith   ierr = ISDestroy(istmp); CHKERRQ(ierr);
6101eb62cbbSBarry Smith 
6111eb62cbbSBarry Smith   /* wait on sends */
6121eb62cbbSBarry Smith   if (nsends) {
6130452661fSBarry Smith     send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));
61478b31e54SBarry Smith     CHKPTRQ(send_status);
6151eb62cbbSBarry Smith     MPI_Waitall(nsends,send_waits,send_status);
6160452661fSBarry Smith     PetscFree(send_status);
6171eb62cbbSBarry Smith   }
6180452661fSBarry Smith   PetscFree(send_waits); PetscFree(svalues);
6191eb62cbbSBarry Smith 
6201eb62cbbSBarry Smith   return 0;
6211eb62cbbSBarry Smith }
6221eb62cbbSBarry Smith 
6235615d1e5SSatish Balay #undef __FUNC__
6245615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIAIJ"
6258f6be9afSLois Curfman McInnes int MatMult_MPIAIJ(Mat A,Vec xx,Vec yy)
6261eb62cbbSBarry Smith {
627416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
628fbd6ef76SBarry Smith   int        ierr,nt;
629416022c9SBarry Smith 
630a2ce50c7SBarry Smith   ierr = VecGetLocalSize(xx,&nt);  CHKERRQ(ierr);
631fbd6ef76SBarry Smith   if (nt != a->n) {
632f40265daSLois Curfman McInnes     SETERRQ(1,0,"Incompatible partition of A and xx");
633fbd6ef76SBarry Smith   }
63443a90d84SBarry Smith   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr);
63538597bd4SSatish Balay   ierr = (*a->A->ops.mult)(a->A,xx,yy); CHKERRQ(ierr);
63643a90d84SBarry Smith   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx); CHKERRQ(ierr);
63738597bd4SSatish Balay   ierr = (*a->B->ops.multadd)(a->B,a->lvec,yy,yy); CHKERRQ(ierr);
6381eb62cbbSBarry Smith   return 0;
6391eb62cbbSBarry Smith }
6401eb62cbbSBarry Smith 
6415615d1e5SSatish Balay #undef __FUNC__
6425615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIAIJ"
6438f6be9afSLois Curfman McInnes int MatMultAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
644da3a660dSBarry Smith {
645416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
646da3a660dSBarry Smith   int        ierr;
64743a90d84SBarry Smith   ierr = VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
64838597bd4SSatish Balay   ierr = (*a->A->ops.multadd)(a->A,xx,yy,zz); CHKERRQ(ierr);
64943a90d84SBarry Smith   ierr = VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
65038597bd4SSatish Balay   ierr = (*a->B->ops.multadd)(a->B,a->lvec,zz,zz); CHKERRQ(ierr);
651da3a660dSBarry Smith   return 0;
652da3a660dSBarry Smith }
653da3a660dSBarry Smith 
6545615d1e5SSatish Balay #undef __FUNC__
6555615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIAIJ"
6568f6be9afSLois Curfman McInnes int MatMultTrans_MPIAIJ(Mat A,Vec xx,Vec yy)
657da3a660dSBarry Smith {
658416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
659da3a660dSBarry Smith   int        ierr;
660da3a660dSBarry Smith 
661da3a660dSBarry Smith   /* do nondiagonal part */
66238597bd4SSatish Balay   ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr);
663da3a660dSBarry Smith   /* send it on its way */
664537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
665da3a660dSBarry Smith   /* do local part */
66638597bd4SSatish Balay   ierr = (*a->A->ops.multtrans)(a->A,xx,yy); CHKERRQ(ierr);
667da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
668da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
669da3a660dSBarry Smith   /* but is not perhaps always true. */
670537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
671da3a660dSBarry Smith   return 0;
672da3a660dSBarry Smith }
673da3a660dSBarry Smith 
6745615d1e5SSatish Balay #undef __FUNC__
6755615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIAIJ"
6768f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIAIJ(Mat A,Vec xx,Vec yy,Vec zz)
677da3a660dSBarry Smith {
678416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
679da3a660dSBarry Smith   int        ierr;
680da3a660dSBarry Smith 
681da3a660dSBarry Smith   /* do nondiagonal part */
68238597bd4SSatish Balay   ierr = (*a->B->ops.multtrans)(a->B,xx,a->lvec); CHKERRQ(ierr);
683da3a660dSBarry Smith   /* send it on its way */
684537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
685da3a660dSBarry Smith   /* do local part */
68638597bd4SSatish Balay   ierr = (*a->A->ops.multtransadd)(a->A,xx,yy,zz); CHKERRQ(ierr);
687da3a660dSBarry Smith   /* receive remote parts: note this assumes the values are not actually */
688da3a660dSBarry Smith   /* inserted in yy until the next line, which is true for my implementation*/
689da3a660dSBarry Smith   /* but is not perhaps always true. */
6900a198c4cSBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr);
691da3a660dSBarry Smith   return 0;
692da3a660dSBarry Smith }
693da3a660dSBarry Smith 
6941eb62cbbSBarry Smith /*
6951eb62cbbSBarry Smith   This only works correctly for square matrices where the subblock A->A is the
6961eb62cbbSBarry Smith    diagonal block
6971eb62cbbSBarry Smith */
6985615d1e5SSatish Balay #undef __FUNC__
6995615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIAIJ"
7008f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIAIJ(Mat A,Vec v)
7011eb62cbbSBarry Smith {
702416022c9SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
70344cd7ae7SLois Curfman McInnes   if (a->M != a->N)
704e3372554SBarry Smith     SETERRQ(1,0,"Supports only square matrix where A->A is diag block");
7055baf8537SBarry Smith   if (a->rstart != a->cstart || a->rend != a->cend) {
706e3372554SBarry Smith     SETERRQ(1,0,"row partition must equal col partition");  }
707416022c9SBarry Smith   return MatGetDiagonal(a->A,v);
7081eb62cbbSBarry Smith }
7091eb62cbbSBarry Smith 
7105615d1e5SSatish Balay #undef __FUNC__
7115615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIAIJ"
7128f6be9afSLois Curfman McInnes int MatScale_MPIAIJ(Scalar *aa,Mat A)
713052efed2SBarry Smith {
714052efed2SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
715052efed2SBarry Smith   int        ierr;
716052efed2SBarry Smith   ierr = MatScale(aa,a->A); CHKERRQ(ierr);
717052efed2SBarry Smith   ierr = MatScale(aa,a->B); CHKERRQ(ierr);
718052efed2SBarry Smith   return 0;
719052efed2SBarry Smith }
720052efed2SBarry Smith 
7215615d1e5SSatish Balay #undef __FUNC__
7225eea60f9SBarry Smith #define __FUNC__ "MatDestroy_MPIAIJ" /* ADIC Ignore */
7238f6be9afSLois Curfman McInnes int MatDestroy_MPIAIJ(PetscObject obj)
7241eb62cbbSBarry Smith {
7251eb62cbbSBarry Smith   Mat        mat = (Mat) obj;
72644a69424SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
7271eb62cbbSBarry Smith   int        ierr;
728a5a9c739SBarry Smith #if defined(PETSC_LOG)
7296e35fa57SLois Curfman McInnes   PLogObjectState(obj,"Rows=%d, Cols=%d",aij->M,aij->N);
730a5a9c739SBarry Smith #endif
7310452661fSBarry Smith   PetscFree(aij->rowners);
73278b31e54SBarry Smith   ierr = MatDestroy(aij->A); CHKERRQ(ierr);
73378b31e54SBarry Smith   ierr = MatDestroy(aij->B); CHKERRQ(ierr);
7340452661fSBarry Smith   if (aij->colmap) PetscFree(aij->colmap);
7350452661fSBarry Smith   if (aij->garray) PetscFree(aij->garray);
7361eb62cbbSBarry Smith   if (aij->lvec)   VecDestroy(aij->lvec);
737a56f8943SBarry Smith   if (aij->Mvctx)  VecScatterDestroy(aij->Mvctx);
7387a0afa10SBarry Smith   if (aij->rowvalues) PetscFree(aij->rowvalues);
7390452661fSBarry Smith   PetscFree(aij);
74090f02eecSBarry Smith   if (mat->mapping) {
74190f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
74290f02eecSBarry Smith   }
743a5a9c739SBarry Smith   PLogObjectDestroy(mat);
7440452661fSBarry Smith   PetscHeaderDestroy(mat);
7451eb62cbbSBarry Smith   return 0;
7461eb62cbbSBarry Smith }
747ee50ffe9SBarry Smith 
7485615d1e5SSatish Balay #undef __FUNC__
7495eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_Binary" /* ADIC Ignore */
7508f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_Binary(Mat mat,Viewer viewer)
7511eb62cbbSBarry Smith {
752416022c9SBarry Smith   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
753416022c9SBarry Smith   int         ierr;
754416022c9SBarry Smith 
75517699dbbSLois Curfman McInnes   if (aij->size == 1) {
756416022c9SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
757416022c9SBarry Smith   }
758e3372554SBarry Smith   else SETERRQ(1,0,"Only uniprocessor output supported");
759416022c9SBarry Smith   return 0;
760416022c9SBarry Smith }
761416022c9SBarry Smith 
7625615d1e5SSatish Balay #undef __FUNC__
7635eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ_ASCIIorDraworMatlab" /* ADIC Ignore */
7648f6be9afSLois Curfman McInnes extern int MatView_MPIAIJ_ASCIIorDraworMatlab(Mat mat,Viewer viewer)
765416022c9SBarry Smith {
76644a69424SLois Curfman McInnes   Mat_MPIAIJ  *aij = (Mat_MPIAIJ *) mat->data;
767dbb450caSBarry Smith   Mat_SeqAIJ* C = (Mat_SeqAIJ*)aij->A->data;
768a56f8943SBarry Smith   int         ierr, format,shift = C->indexshift,rank;
769d636dbe3SBarry Smith   FILE        *fd;
77019bcc07fSBarry Smith   ViewerType  vtype;
771416022c9SBarry Smith 
77219bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
77319bcc07fSBarry Smith   if (vtype  == ASCII_FILES_VIEWER || vtype == ASCII_FILE_VIEWER) {
77490ace30eSBarry Smith     ierr = ViewerGetFormat(viewer,&format);
7750a198c4cSBarry Smith     if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
7764e220ebcSLois Curfman McInnes       MatInfo info;
7774e220ebcSLois Curfman McInnes       int     flg;
778a56f8943SBarry Smith       MPI_Comm_rank(mat->comm,&rank);
77990ace30eSBarry Smith       ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
7804e220ebcSLois Curfman McInnes       ierr = MatGetInfo(mat,MAT_LOCAL,&info);
78195e01e2fSLois Curfman McInnes       ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr);
78277c4ece6SBarry Smith       PetscSequentialPhaseBegin(mat->comm,1);
78395e01e2fSLois Curfman McInnes       if (flg) fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, not using I-node routines\n",
7844e220ebcSLois Curfman McInnes          rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
78595e01e2fSLois Curfman McInnes       else fprintf(fd,"[%d] Local rows %d nz %d nz alloced %d mem %d, using I-node routines\n",
7864e220ebcSLois Curfman McInnes          rank,aij->m,(int)info.nz_used,(int)info.nz_allocated,(int)info.memory);
7874e220ebcSLois Curfman McInnes       ierr = MatGetInfo(aij->A,MAT_LOCAL,&info);
7884e220ebcSLois Curfman McInnes       fprintf(fd,"[%d] on-diagonal part: nz %d \n",rank,(int)info.nz_used);
7894e220ebcSLois Curfman McInnes       ierr = MatGetInfo(aij->B,MAT_LOCAL,&info);
7904e220ebcSLois Curfman McInnes       fprintf(fd,"[%d] off-diagonal part: nz %d \n",rank,(int)info.nz_used);
791a56f8943SBarry Smith       fflush(fd);
79277c4ece6SBarry Smith       PetscSequentialPhaseEnd(mat->comm,1);
793a40aa06bSLois Curfman McInnes       ierr = VecScatterView(aij->Mvctx,viewer); CHKERRQ(ierr);
794416022c9SBarry Smith       return 0;
795416022c9SBarry Smith     }
7960a198c4cSBarry Smith     else if (format == VIEWER_FORMAT_ASCII_INFO) {
79708480c60SBarry Smith       return 0;
79808480c60SBarry Smith     }
799416022c9SBarry Smith   }
800416022c9SBarry Smith 
80119bcc07fSBarry Smith   if (vtype == DRAW_VIEWER) {
80219bcc07fSBarry Smith     Draw       draw;
80319bcc07fSBarry Smith     PetscTruth isnull;
80419bcc07fSBarry Smith     ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
80519bcc07fSBarry Smith     ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
80619bcc07fSBarry Smith   }
80719bcc07fSBarry Smith 
80819bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER) {
80990ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
81077c4ece6SBarry Smith     PetscSequentialPhaseBegin(mat->comm,1);
811d13ab20cSBarry Smith     fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n",
81217699dbbSLois Curfman McInnes            aij->rank,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart,
8131eb62cbbSBarry Smith            aij->cend);
81478b31e54SBarry Smith     ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
81578b31e54SBarry Smith     ierr = MatView(aij->B,viewer); CHKERRQ(ierr);
816d13ab20cSBarry Smith     fflush(fd);
81777c4ece6SBarry Smith     PetscSequentialPhaseEnd(mat->comm,1);
818d13ab20cSBarry Smith   }
819416022c9SBarry Smith   else {
820a56f8943SBarry Smith     int size = aij->size;
821a56f8943SBarry Smith     rank = aij->rank;
82217699dbbSLois Curfman McInnes     if (size == 1) {
82378b31e54SBarry Smith       ierr = MatView(aij->A,viewer); CHKERRQ(ierr);
82495373324SBarry Smith     }
82595373324SBarry Smith     else {
82695373324SBarry Smith       /* assemble the entire matrix onto first processor. */
82795373324SBarry Smith       Mat         A;
828ec8511deSBarry Smith       Mat_SeqAIJ *Aloc;
8292eb8c8abSBarry Smith       int         M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct;
83095373324SBarry Smith       Scalar      *a;
8312ee70a88SLois Curfman McInnes 
83217699dbbSLois Curfman McInnes       if (!rank) {
833b4fd4287SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);
834c451ab25SLois Curfman McInnes                CHKERRQ(ierr);
83595373324SBarry Smith       }
83695373324SBarry Smith       else {
837b4fd4287SBarry Smith         ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,PETSC_NULL,0,PETSC_NULL,&A);
838c451ab25SLois Curfman McInnes                CHKERRQ(ierr);
83995373324SBarry Smith       }
840464493b3SBarry Smith       PLogObjectParent(mat,A);
841416022c9SBarry Smith 
84295373324SBarry Smith       /* copy over the A part */
843ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->A->data;
8442ee70a88SLois Curfman McInnes       m = Aloc->m; ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
84595373324SBarry Smith       row = aij->rstart;
846dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += aij->cstart + shift;}
84795373324SBarry Smith       for ( i=0; i<m; i++ ) {
848416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,INSERT_VALUES);CHKERRQ(ierr);
84995373324SBarry Smith         row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
85095373324SBarry Smith       }
8512ee70a88SLois Curfman McInnes       aj = Aloc->j;
852dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= aij->cstart + shift;}
85395373324SBarry Smith 
85495373324SBarry Smith       /* copy over the B part */
855ec8511deSBarry Smith       Aloc = (Mat_SeqAIJ*) aij->B->data;
8562ee70a88SLois Curfman McInnes       m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
85795373324SBarry Smith       row = aij->rstart;
8580452661fSBarry Smith       ct = cols = (int *) PetscMalloc( (ai[m]+1)*sizeof(int) ); CHKPTRQ(cols);
859dbb450caSBarry Smith       for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = aij->garray[aj[i]+shift];}
86095373324SBarry Smith       for ( i=0; i<m; i++ ) {
861416022c9SBarry Smith         ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,INSERT_VALUES);CHKERRQ(ierr);
86295373324SBarry Smith         row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
86395373324SBarry Smith       }
8640452661fSBarry Smith       PetscFree(ct);
8656d4a8577SBarry Smith       ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8666d4a8577SBarry Smith       ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
86717699dbbSLois Curfman McInnes       if (!rank) {
86878b31e54SBarry Smith         ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERRQ(ierr);
86995373324SBarry Smith       }
87078b31e54SBarry Smith       ierr = MatDestroy(A); CHKERRQ(ierr);
87195373324SBarry Smith     }
87295373324SBarry Smith   }
8731eb62cbbSBarry Smith   return 0;
8741eb62cbbSBarry Smith }
8751eb62cbbSBarry Smith 
8765615d1e5SSatish Balay #undef __FUNC__
8775eea60f9SBarry Smith #define __FUNC__ "MatView_MPIAIJ" /* ADIC Ignore */
8788f6be9afSLois Curfman McInnes int MatView_MPIAIJ(PetscObject obj,Viewer viewer)
879416022c9SBarry Smith {
880416022c9SBarry Smith   Mat         mat = (Mat) obj;
881416022c9SBarry Smith   int         ierr;
88219bcc07fSBarry Smith   ViewerType  vtype;
883416022c9SBarry Smith 
88419bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
88519bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER ||
88619bcc07fSBarry Smith       vtype == DRAW_VIEWER       || vtype == MATLAB_VIEWER) {
887d7e8b826SBarry Smith     ierr = MatView_MPIAIJ_ASCIIorDraworMatlab(mat,viewer); CHKERRQ(ierr);
888416022c9SBarry Smith   }
88919bcc07fSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
890416022c9SBarry Smith     return MatView_MPIAIJ_Binary(mat,viewer);
891416022c9SBarry Smith   }
892416022c9SBarry Smith   return 0;
893416022c9SBarry Smith }
894416022c9SBarry Smith 
8951eb62cbbSBarry Smith /*
8961eb62cbbSBarry Smith     This has to provide several versions.
8971eb62cbbSBarry Smith 
8981eb62cbbSBarry Smith      2) a) use only local smoothing updating outer values only once.
8991eb62cbbSBarry Smith         b) local smoothing updating outer values each inner iteration
900d6dfbf8fSBarry Smith      3) color updating out values betwen colors.
9011eb62cbbSBarry Smith */
9025615d1e5SSatish Balay #undef __FUNC__
9035615d1e5SSatish Balay #define __FUNC__ "MatRelax_MPIAIJ"
9048f6be9afSLois Curfman McInnes int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,MatSORType flag,
905dbb450caSBarry Smith                            double fshift,int its,Vec xx)
9068a729477SBarry Smith {
90744a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
908d6dfbf8fSBarry Smith   Mat        AA = mat->A, BB = mat->B;
909ec8511deSBarry Smith   Mat_SeqAIJ *A = (Mat_SeqAIJ *) AA->data, *B = (Mat_SeqAIJ *)BB->data;
910c16cb8f2SBarry Smith   Scalar     *b,*x,*xs,*ls,d,*v,sum;
9116abc6512SBarry Smith   int        ierr,*idx, *diag;
912416022c9SBarry Smith   int        n = mat->n, m = mat->m, i,shift = A->indexshift;
9138a729477SBarry Smith 
914d6dfbf8fSBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls);
915dbb450caSBarry Smith   xs = x + shift; /* shift by one for index start of 1 */
916dbb450caSBarry Smith   ls = ls + shift;
917ec8511deSBarry Smith   if (!A->diag) {if ((ierr = MatMarkDiag_SeqAIJ(AA))) return ierr;}
918d6dfbf8fSBarry Smith   diag = A->diag;
919c16cb8f2SBarry Smith   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
920da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
92138597bd4SSatish Balay       return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);
922da3a660dSBarry Smith     }
92343a90d84SBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
92478b31e54SBarry Smith     CHKERRQ(ierr);
92543a90d84SBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
92678b31e54SBarry Smith     CHKERRQ(ierr);
927d6dfbf8fSBarry Smith     while (its--) {
928d6dfbf8fSBarry Smith       /* go down through the rows */
929d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
930d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
931dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
932dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
933d6dfbf8fSBarry Smith         sum  = b[i];
934d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
935dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
936d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
937dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
938dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
939d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
94055a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
941d6dfbf8fSBarry Smith       }
942d6dfbf8fSBarry Smith       /* come up through the rows */
943d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
944d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
945dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
946dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
947d6dfbf8fSBarry Smith         sum  = b[i];
948d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
949dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
950d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
951dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
952dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
953d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
95455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
955d6dfbf8fSBarry Smith       }
956d6dfbf8fSBarry Smith     }
957d6dfbf8fSBarry Smith   }
958d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_FORWARD_SWEEP){
959da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
96038597bd4SSatish Balay       return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);
961da3a660dSBarry Smith     }
96243a90d84SBarry Smith     ierr=VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
96378b31e54SBarry Smith     CHKERRQ(ierr);
96443a90d84SBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
96578b31e54SBarry Smith     CHKERRQ(ierr);
966d6dfbf8fSBarry Smith     while (its--) {
967d6dfbf8fSBarry Smith       for ( i=0; i<m; i++ ) {
968d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
969dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
970dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
971d6dfbf8fSBarry Smith         sum  = b[i];
972d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
973dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
974d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
975dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
976dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
977d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
97855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
979d6dfbf8fSBarry Smith       }
980d6dfbf8fSBarry Smith     }
981d6dfbf8fSBarry Smith   }
982d6dfbf8fSBarry Smith   else if (flag & SOR_LOCAL_BACKWARD_SWEEP){
983da3a660dSBarry Smith     if (flag & SOR_ZERO_INITIAL_GUESS) {
98438597bd4SSatish Balay       return (*mat->A->ops.relax)(mat->A,bb,omega,flag,fshift,its,xx);
985da3a660dSBarry Smith     }
98643a90d84SBarry Smith     ierr = VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,
98778b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
98843a90d84SBarry Smith     ierr = VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,
98978b31e54SBarry Smith                             mat->Mvctx); CHKERRQ(ierr);
990d6dfbf8fSBarry Smith     while (its--) {
991d6dfbf8fSBarry Smith       for ( i=m-1; i>-1; i-- ) {
992d6dfbf8fSBarry Smith         n    = A->i[i+1] - A->i[i];
993dbb450caSBarry Smith         idx  = A->j + A->i[i] + shift;
994dbb450caSBarry Smith         v    = A->a + A->i[i] + shift;
995d6dfbf8fSBarry Smith         sum  = b[i];
996d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
997dbb450caSBarry Smith         d    = fshift + A->a[diag[i]+shift];
998d6dfbf8fSBarry Smith         n    = B->i[i+1] - B->i[i];
999dbb450caSBarry Smith         idx  = B->j + B->i[i] + shift;
1000dbb450caSBarry Smith         v    = B->a + B->i[i] + shift;
1001d6dfbf8fSBarry Smith         SPARSEDENSEMDOT(sum,ls,v,idx,n);
100255a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + A->a[diag[i]+shift]*x[i])/d;
1003d6dfbf8fSBarry Smith       }
1004d6dfbf8fSBarry Smith     }
1005d6dfbf8fSBarry Smith   }
1006c16cb8f2SBarry Smith   else {
1007e3372554SBarry Smith     SETERRQ(1,0,"Option not supported");
1008c16cb8f2SBarry Smith   }
10098a729477SBarry Smith   return 0;
10108a729477SBarry Smith }
1011a66be287SLois Curfman McInnes 
10125615d1e5SSatish Balay #undef __FUNC__
10135eea60f9SBarry Smith #define __FUNC__ "MatGetInfo_MPIAIJ" /* ADIC Ignore */
10148f6be9afSLois Curfman McInnes int MatGetInfo_MPIAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1015a66be287SLois Curfman McInnes {
1016a66be287SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1017a66be287SLois Curfman McInnes   Mat        A = mat->A, B = mat->B;
10184e220ebcSLois Curfman McInnes   int        ierr;
10194e220ebcSLois Curfman McInnes   double     isend[5], irecv[5];
1020a66be287SLois Curfman McInnes 
10214e220ebcSLois Curfman McInnes   info->rows_global    = (double)mat->M;
10224e220ebcSLois Curfman McInnes   info->columns_global = (double)mat->N;
10234e220ebcSLois Curfman McInnes   info->rows_local     = (double)mat->m;
10244e220ebcSLois Curfman McInnes   info->columns_local  = (double)mat->N;
10254e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10264e220ebcSLois Curfman McInnes   ierr = MatGetInfo(A,MAT_LOCAL,info); CHKERRQ(ierr);
10274e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
10284e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
10294e220ebcSLois Curfman McInnes   ierr = MatGetInfo(B,MAT_LOCAL,info); CHKERRQ(ierr);
10304e220ebcSLois Curfman McInnes   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
10314e220ebcSLois Curfman McInnes   isend[3] += info->memory;  isend[4] += info->mallocs;
1032a66be287SLois Curfman McInnes   if (flag == MAT_LOCAL) {
10334e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
10344e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
10354e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
10364e220ebcSLois Curfman McInnes     info->memory       = isend[3];
10374e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
1038a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
10394e220ebcSLois Curfman McInnes     MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,matin->comm);
10404e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
10414e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
10424e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
10434e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
10444e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
1045a66be287SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
10464e220ebcSLois Curfman McInnes     MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,matin->comm);
10474e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
10484e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
10494e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
10504e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
10514e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
1052a66be287SLois Curfman McInnes   }
10534e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
10544e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
10554e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
10564e220ebcSLois Curfman McInnes 
1057a66be287SLois Curfman McInnes   return 0;
1058a66be287SLois Curfman McInnes }
1059a66be287SLois Curfman McInnes 
1060299609e3SLois Curfman McInnes extern int MatLUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,Mat*);
1061299609e3SLois Curfman McInnes extern int MatLUFactorNumeric_MPIAIJ(Mat,Mat*);
1062299609e3SLois Curfman McInnes extern int MatLUFactor_MPIAIJ(Mat,IS,IS,double);
1063299609e3SLois Curfman McInnes extern int MatILUFactorSymbolic_MPIAIJ(Mat,IS,IS,double,int,Mat *);
1064299609e3SLois Curfman McInnes extern int MatSolve_MPIAIJ(Mat,Vec,Vec);
1065299609e3SLois Curfman McInnes extern int MatSolveAdd_MPIAIJ(Mat,Vec,Vec,Vec);
1066299609e3SLois Curfman McInnes extern int MatSolveTrans_MPIAIJ(Mat,Vec,Vec);
1067299609e3SLois Curfman McInnes extern int MatSolveTransAdd_MPIAIJ(Mat,Vec,Vec,Vec);
1068299609e3SLois Curfman McInnes 
10695615d1e5SSatish Balay #undef __FUNC__
10705eea60f9SBarry Smith #define __FUNC__ "MatSetOption_MPIAIJ" /* ADIC Ignore */
10718f6be9afSLois Curfman McInnes int MatSetOption_MPIAIJ(Mat A,MatOption op)
1072c74985f6SBarry Smith {
1073c0bbcb79SLois Curfman McInnes   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1074c74985f6SBarry Smith 
10756d4a8577SBarry Smith   if (op == MAT_NO_NEW_NONZERO_LOCATIONS ||
10766d4a8577SBarry Smith       op == MAT_YES_NEW_NONZERO_LOCATIONS ||
10776da5968aSLois Curfman McInnes       op == MAT_COLUMNS_UNSORTED ||
1078c2653b3dSLois Curfman McInnes       op == MAT_COLUMNS_SORTED ||
107996854ed6SLois Curfman McInnes       op == MAT_NEW_NONZERO_ALLOCATION_ERROR ||
1080c2653b3dSLois Curfman McInnes       op == MAT_NEW_NONZERO_LOCATION_ERROR) {
1081b1fbbac0SLois Curfman McInnes         MatSetOption(a->A,op);
1082b1fbbac0SLois Curfman McInnes         MatSetOption(a->B,op);
1083b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROW_ORIENTED) {
1084aeafbbfcSLois Curfman McInnes     a->roworiented = 1;
1085c0bbcb79SLois Curfman McInnes     MatSetOption(a->A,op);
1086c0bbcb79SLois Curfman McInnes     MatSetOption(a->B,op);
1087b1fbbac0SLois Curfman McInnes   } else if (op == MAT_ROWS_SORTED ||
10886da5968aSLois Curfman McInnes              op == MAT_ROWS_UNSORTED ||
10896d4a8577SBarry Smith              op == MAT_SYMMETRIC ||
10906d4a8577SBarry Smith              op == MAT_STRUCTURALLY_SYMMETRIC ||
10916d4a8577SBarry Smith              op == MAT_YES_NEW_DIAGONALS)
109294a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_MPIAIJ:Option ignored\n");
10936d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED) {
10944b0e389bSBarry Smith     a->roworiented = 0;
10954b0e389bSBarry Smith     MatSetOption(a->A,op);
10964b0e389bSBarry Smith     MatSetOption(a->B,op);
10972b362799SSatish Balay   } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) {
109890f02eecSBarry Smith     a->donotstash = 1;
109990f02eecSBarry Smith   } else if (op == MAT_NO_NEW_DIAGONALS)
1100e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");}
1101c0bbcb79SLois Curfman McInnes   else
1102e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
1103c74985f6SBarry Smith   return 0;
1104c74985f6SBarry Smith }
1105c74985f6SBarry Smith 
11065615d1e5SSatish Balay #undef __FUNC__
11075eea60f9SBarry Smith #define __FUNC__ "MatGetSize_MPIAIJ" /* ADIC Ignore */
11088f6be9afSLois Curfman McInnes int MatGetSize_MPIAIJ(Mat matin,int *m,int *n)
1109c74985f6SBarry Smith {
111044a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1111c74985f6SBarry Smith   *m = mat->M; *n = mat->N;
1112c74985f6SBarry Smith   return 0;
1113c74985f6SBarry Smith }
1114c74985f6SBarry Smith 
11155615d1e5SSatish Balay #undef __FUNC__
11165eea60f9SBarry Smith #define __FUNC__ "MatGetLocalSize_MPIAIJ" /* ADIC Ignore */
11178f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIAIJ(Mat matin,int *m,int *n)
1118c74985f6SBarry Smith {
111944a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1120b7c46309SBarry Smith   *m = mat->m; *n = mat->N;
1121c74985f6SBarry Smith   return 0;
1122c74985f6SBarry Smith }
1123c74985f6SBarry Smith 
11245615d1e5SSatish Balay #undef __FUNC__
11255eea60f9SBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAIJ" /* ADIC Ignore */
11268f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIAIJ(Mat matin,int *m,int *n)
1127c74985f6SBarry Smith {
112844a69424SLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
1129c74985f6SBarry Smith   *m = mat->rstart; *n = mat->rend;
1130c74985f6SBarry Smith   return 0;
1131c74985f6SBarry Smith }
1132c74985f6SBarry Smith 
11336d84be18SBarry Smith extern int MatGetRow_SeqAIJ(Mat,int,int*,int**,Scalar**);
11346d84be18SBarry Smith extern int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,Scalar**);
11356d84be18SBarry Smith 
11365615d1e5SSatish Balay #undef __FUNC__
11375615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIAIJ"
11386d84be18SBarry Smith int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v)
113939e00950SLois Curfman McInnes {
1140154123eaSLois Curfman McInnes   Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data;
114170f0671dSBarry Smith   Scalar     *vworkA, *vworkB, **pvA, **pvB,*v_p;
1142154123eaSLois Curfman McInnes   int        i, ierr, *cworkA, *cworkB, **pcA, **pcB, cstart = mat->cstart;
1143154123eaSLois Curfman McInnes   int        nztot, nzA, nzB, lrow, rstart = mat->rstart, rend = mat->rend;
114470f0671dSBarry Smith   int        *cmap, *idx_p;
114539e00950SLois Curfman McInnes 
1146e3372554SBarry Smith   if (mat->getrowactive == PETSC_TRUE) SETERRQ(1,0,"Already active");
11477a0afa10SBarry Smith   mat->getrowactive = PETSC_TRUE;
11487a0afa10SBarry Smith 
114970f0671dSBarry Smith   if (!mat->rowvalues && (idx || v)) {
11507a0afa10SBarry Smith     /*
11517a0afa10SBarry Smith         allocate enough space to hold information from the longest row.
11527a0afa10SBarry Smith     */
11537a0afa10SBarry Smith     Mat_SeqAIJ *Aa = (Mat_SeqAIJ *) mat->A->data,*Ba = (Mat_SeqAIJ *) mat->B->data;
1154c16cb8f2SBarry Smith     int     max = 1,m = mat->m,tmp;
1155c16cb8f2SBarry Smith     for ( i=0; i<m; i++ ) {
11567a0afa10SBarry Smith       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i];
11577a0afa10SBarry Smith       if (max < tmp) { max = tmp; }
11587a0afa10SBarry Smith     }
11597a0afa10SBarry Smith     mat->rowvalues = (Scalar *) PetscMalloc( max*(sizeof(int)+sizeof(Scalar)));
11607a0afa10SBarry Smith     CHKPTRQ(mat->rowvalues);
11617a0afa10SBarry Smith     mat->rowindices = (int *) (mat->rowvalues + max);
11627a0afa10SBarry Smith   }
11637a0afa10SBarry Smith 
1164e3372554SBarry Smith   if (row < rstart || row >= rend) SETERRQ(1,0,"Only local rows")
1165abc0e9e4SLois Curfman McInnes   lrow = row - rstart;
116639e00950SLois Curfman McInnes 
1167154123eaSLois Curfman McInnes   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1168154123eaSLois Curfman McInnes   if (!v)   {pvA = 0; pvB = 0;}
1169154123eaSLois Curfman McInnes   if (!idx) {pcA = 0; if (!v) pcB = 0;}
117038597bd4SSatish Balay   ierr = (*mat->A->ops.getrow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
117138597bd4SSatish Balay   ierr = (*mat->B->ops.getrow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
1172154123eaSLois Curfman McInnes   nztot = nzA + nzB;
1173154123eaSLois Curfman McInnes 
117470f0671dSBarry Smith   cmap  = mat->garray;
1175154123eaSLois Curfman McInnes   if (v  || idx) {
1176154123eaSLois Curfman McInnes     if (nztot) {
1177154123eaSLois Curfman McInnes       /* Sort by increasing column numbers, assuming A and B already sorted */
117870f0671dSBarry Smith       int imark = -1;
1179154123eaSLois Curfman McInnes       if (v) {
118070f0671dSBarry Smith         *v = v_p = mat->rowvalues;
118139e00950SLois Curfman McInnes         for ( i=0; i<nzB; i++ ) {
118270f0671dSBarry Smith           if (cmap[cworkB[i]] < cstart)   v_p[i] = vworkB[i];
1183154123eaSLois Curfman McInnes           else break;
1184154123eaSLois Curfman McInnes         }
1185154123eaSLois Curfman McInnes         imark = i;
118670f0671dSBarry Smith         for ( i=0; i<nzA; i++ )     v_p[imark+i] = vworkA[i];
118770f0671dSBarry Smith         for ( i=imark; i<nzB; i++ ) v_p[nzA+i]   = vworkB[i];
1188154123eaSLois Curfman McInnes       }
1189154123eaSLois Curfman McInnes       if (idx) {
119070f0671dSBarry Smith         *idx = idx_p = mat->rowindices;
119170f0671dSBarry Smith         if (imark > -1) {
119270f0671dSBarry Smith           for ( i=0; i<imark; i++ ) {
119370f0671dSBarry Smith             idx_p[i] = cmap[cworkB[i]];
119470f0671dSBarry Smith           }
119570f0671dSBarry Smith         } else {
1196154123eaSLois Curfman McInnes           for ( i=0; i<nzB; i++ ) {
119770f0671dSBarry Smith             if (cmap[cworkB[i]] < cstart)   idx_p[i] = cmap[cworkB[i]];
1198154123eaSLois Curfman McInnes             else break;
1199154123eaSLois Curfman McInnes           }
1200154123eaSLois Curfman McInnes           imark = i;
120170f0671dSBarry Smith         }
120270f0671dSBarry Smith         for ( i=0; i<nzA; i++ )     idx_p[imark+i] = cstart + cworkA[i];
120370f0671dSBarry Smith         for ( i=imark; i<nzB; i++ ) idx_p[nzA+i]   = cmap[cworkB[i]];
120439e00950SLois Curfman McInnes       }
120539e00950SLois Curfman McInnes     }
12061ca473b0SSatish Balay     else {
12071ca473b0SSatish Balay       if (idx) *idx = 0;
12081ca473b0SSatish Balay       if (v)   *v   = 0;
12091ca473b0SSatish Balay     }
1210154123eaSLois Curfman McInnes   }
121139e00950SLois Curfman McInnes   *nz = nztot;
121238597bd4SSatish Balay   ierr = (*mat->A->ops.restorerow)(mat->A,lrow,&nzA,pcA,pvA); CHKERRQ(ierr);
121338597bd4SSatish Balay   ierr = (*mat->B->ops.restorerow)(mat->B,lrow,&nzB,pcB,pvB); CHKERRQ(ierr);
121439e00950SLois Curfman McInnes   return 0;
121539e00950SLois Curfman McInnes }
121639e00950SLois Curfman McInnes 
12175615d1e5SSatish Balay #undef __FUNC__
12185eea60f9SBarry Smith #define __FUNC__ "MatRestoreRow_MPIAIJ" /* ADIC Ignore */
12196d84be18SBarry Smith int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v)
122039e00950SLois Curfman McInnes {
12217a0afa10SBarry Smith   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
12227a0afa10SBarry Smith   if (aij->getrowactive == PETSC_FALSE) {
1223e3372554SBarry Smith     SETERRQ(1,0,"MatGetRow not called");
12247a0afa10SBarry Smith   }
12257a0afa10SBarry Smith   aij->getrowactive = PETSC_FALSE;
122639e00950SLois Curfman McInnes   return 0;
122739e00950SLois Curfman McInnes }
122839e00950SLois Curfman McInnes 
12295615d1e5SSatish Balay #undef __FUNC__
12305615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIAIJ"
12318f6be9afSLois Curfman McInnes int MatNorm_MPIAIJ(Mat mat,NormType type,double *norm)
1232855ac2c5SLois Curfman McInnes {
1233855ac2c5SLois Curfman McInnes   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
1234ec8511deSBarry Smith   Mat_SeqAIJ *amat = (Mat_SeqAIJ*) aij->A->data, *bmat = (Mat_SeqAIJ*) aij->B->data;
1235416022c9SBarry Smith   int        ierr, i, j, cstart = aij->cstart,shift = amat->indexshift;
1236416022c9SBarry Smith   double     sum = 0.0;
123704ca555eSLois Curfman McInnes   Scalar     *v;
123804ca555eSLois Curfman McInnes 
123917699dbbSLois Curfman McInnes   if (aij->size == 1) {
124014183eadSLois Curfman McInnes     ierr =  MatNorm(aij->A,type,norm); CHKERRQ(ierr);
124137fa93a5SLois Curfman McInnes   } else {
124204ca555eSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
124304ca555eSLois Curfman McInnes       v = amat->a;
124404ca555eSLois Curfman McInnes       for (i=0; i<amat->nz; i++ ) {
124504ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
124604ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
124704ca555eSLois Curfman McInnes #else
124804ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
124904ca555eSLois Curfman McInnes #endif
125004ca555eSLois Curfman McInnes       }
125104ca555eSLois Curfman McInnes       v = bmat->a;
125204ca555eSLois Curfman McInnes       for (i=0; i<bmat->nz; i++ ) {
125304ca555eSLois Curfman McInnes #if defined(PETSC_COMPLEX)
125404ca555eSLois Curfman McInnes         sum += real(conj(*v)*(*v)); v++;
125504ca555eSLois Curfman McInnes #else
125604ca555eSLois Curfman McInnes         sum += (*v)*(*v); v++;
125704ca555eSLois Curfman McInnes #endif
125804ca555eSLois Curfman McInnes       }
12596d84be18SBarry Smith       MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,mat->comm);
126004ca555eSLois Curfman McInnes       *norm = sqrt(*norm);
126104ca555eSLois Curfman McInnes     }
126204ca555eSLois Curfman McInnes     else if (type == NORM_1) { /* max column norm */
126304ca555eSLois Curfman McInnes       double *tmp, *tmp2;
126404ca555eSLois Curfman McInnes       int    *jj, *garray = aij->garray;
1265758f045eSSatish Balay       tmp  = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp);
1266758f045eSSatish Balay       tmp2 = (double *) PetscMalloc( (aij->N+1)*sizeof(double) ); CHKPTRQ(tmp2);
1267cddf8d76SBarry Smith       PetscMemzero(tmp,aij->N*sizeof(double));
126804ca555eSLois Curfman McInnes       *norm = 0.0;
126904ca555eSLois Curfman McInnes       v = amat->a; jj = amat->j;
127004ca555eSLois Curfman McInnes       for ( j=0; j<amat->nz; j++ ) {
1271579c6b6fSBarry Smith         tmp[cstart + *jj++ + shift] += PetscAbsScalar(*v);  v++;
127204ca555eSLois Curfman McInnes       }
127304ca555eSLois Curfman McInnes       v = bmat->a; jj = bmat->j;
127404ca555eSLois Curfman McInnes       for ( j=0; j<bmat->nz; j++ ) {
1275579c6b6fSBarry Smith         tmp[garray[*jj++ + shift]] += PetscAbsScalar(*v); v++;
127604ca555eSLois Curfman McInnes       }
12776d84be18SBarry Smith       MPI_Allreduce(tmp,tmp2,aij->N,MPI_DOUBLE,MPI_SUM,mat->comm);
127804ca555eSLois Curfman McInnes       for ( j=0; j<aij->N; j++ ) {
127904ca555eSLois Curfman McInnes         if (tmp2[j] > *norm) *norm = tmp2[j];
128004ca555eSLois Curfman McInnes       }
12810452661fSBarry Smith       PetscFree(tmp); PetscFree(tmp2);
128204ca555eSLois Curfman McInnes     }
128304ca555eSLois Curfman McInnes     else if (type == NORM_INFINITY) { /* max row norm */
1284515d9167SLois Curfman McInnes       double ntemp = 0.0;
128504ca555eSLois Curfman McInnes       for ( j=0; j<amat->m; j++ ) {
1286dbb450caSBarry Smith         v = amat->a + amat->i[j] + shift;
128704ca555eSLois Curfman McInnes         sum = 0.0;
128804ca555eSLois Curfman McInnes         for ( i=0; i<amat->i[j+1]-amat->i[j]; i++ ) {
1289cddf8d76SBarry Smith           sum += PetscAbsScalar(*v); v++;
129004ca555eSLois Curfman McInnes         }
1291dbb450caSBarry Smith         v = bmat->a + bmat->i[j] + shift;
129204ca555eSLois Curfman McInnes         for ( i=0; i<bmat->i[j+1]-bmat->i[j]; i++ ) {
1293cddf8d76SBarry Smith           sum += PetscAbsScalar(*v); v++;
129404ca555eSLois Curfman McInnes         }
1295515d9167SLois Curfman McInnes         if (sum > ntemp) ntemp = sum;
129604ca555eSLois Curfman McInnes       }
12976d84be18SBarry Smith       MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,mat->comm);
129804ca555eSLois Curfman McInnes     }
129904ca555eSLois Curfman McInnes     else {
1300e3372554SBarry Smith       SETERRQ(1,0,"No support for two norm");
130104ca555eSLois Curfman McInnes     }
130237fa93a5SLois Curfman McInnes   }
1303855ac2c5SLois Curfman McInnes   return 0;
1304855ac2c5SLois Curfman McInnes }
1305855ac2c5SLois Curfman McInnes 
13065615d1e5SSatish Balay #undef __FUNC__
13075615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIAIJ"
13088f6be9afSLois Curfman McInnes int MatTranspose_MPIAIJ(Mat A,Mat *matout)
1309b7c46309SBarry Smith {
1310b7c46309SBarry Smith   Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
1311dbb450caSBarry Smith   Mat_SeqAIJ *Aloc = (Mat_SeqAIJ *) a->A->data;
1312416022c9SBarry Smith   int        ierr,shift = Aloc->indexshift;
1313416022c9SBarry Smith   Mat        B;
1314b7c46309SBarry Smith   int        M = a->M, N = a->N,m,*ai,*aj,row,*cols,i,*ct;
1315b7c46309SBarry Smith   Scalar     *array;
1316b7c46309SBarry Smith 
13173638b69dSLois Curfman McInnes   if (matout == PETSC_NULL && M != N)
1318e3372554SBarry Smith     SETERRQ(1,0,"Square matrix only for in-place");
1319b4fd4287SBarry Smith   ierr = MatCreateMPIAIJ(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,0,PETSC_NULL,0,
1320b4fd4287SBarry Smith          PETSC_NULL,&B); CHKERRQ(ierr);
1321b7c46309SBarry Smith 
1322b7c46309SBarry Smith   /* copy over the A part */
1323ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->A->data;
1324b7c46309SBarry Smith   m = Aloc->m; ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1325b7c46309SBarry Smith   row = a->rstart;
1326dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {aj[i] += a->cstart + shift;}
1327b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1328416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],aj,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1329b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i];
1330b7c46309SBarry Smith   }
1331b7c46309SBarry Smith   aj = Aloc->j;
13324af08d9eSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {aj[i] -= a->cstart + shift;}
1333b7c46309SBarry Smith 
1334b7c46309SBarry Smith   /* copy over the B part */
1335ec8511deSBarry Smith   Aloc = (Mat_SeqAIJ*) a->B->data;
1336b7c46309SBarry Smith   m = Aloc->m;  ai = Aloc->i; aj = Aloc->j; array = Aloc->a;
1337b7c46309SBarry Smith   row = a->rstart;
13380452661fSBarry Smith   ct = cols = (int *) PetscMalloc( (1+ai[m]-shift)*sizeof(int) ); CHKPTRQ(cols);
1339dbb450caSBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) {cols[i] = a->garray[aj[i]+shift];}
1340b7c46309SBarry Smith   for ( i=0; i<m; i++ ) {
1341416022c9SBarry Smith     ierr = MatSetValues(B,ai[i+1]-ai[i],cols,1,&row,array,INSERT_VALUES);CHKERRQ(ierr);
1342b7c46309SBarry Smith     row++; array += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i];
1343b7c46309SBarry Smith   }
13440452661fSBarry Smith   PetscFree(ct);
13456d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13466d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13473638b69dSLois Curfman McInnes   if (matout != PETSC_NULL) {
13480de55854SLois Curfman McInnes     *matout = B;
13490de55854SLois Curfman McInnes   } else {
13500de55854SLois Curfman McInnes     /* This isn't really an in-place transpose .... but free data structures from a */
13510452661fSBarry Smith     PetscFree(a->rowners);
13520de55854SLois Curfman McInnes     ierr = MatDestroy(a->A); CHKERRQ(ierr);
13530de55854SLois Curfman McInnes     ierr = MatDestroy(a->B); CHKERRQ(ierr);
13540452661fSBarry Smith     if (a->colmap) PetscFree(a->colmap);
13550452661fSBarry Smith     if (a->garray) PetscFree(a->garray);
13560de55854SLois Curfman McInnes     if (a->lvec) VecDestroy(a->lvec);
1357a56f8943SBarry Smith     if (a->Mvctx) VecScatterDestroy(a->Mvctx);
13580452661fSBarry Smith     PetscFree(a);
1359416022c9SBarry Smith     PetscMemcpy(A,B,sizeof(struct _Mat));
13600452661fSBarry Smith     PetscHeaderDestroy(B);
13610de55854SLois Curfman McInnes   }
1362b7c46309SBarry Smith   return 0;
1363b7c46309SBarry Smith }
1364b7c46309SBarry Smith 
13655615d1e5SSatish Balay #undef __FUNC__
13665615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_MPIAIJ"
13674b967eb1SSatish Balay int MatDiagonalScale_MPIAIJ(Mat mat,Vec ll,Vec rr)
1368a008b906SSatish Balay {
13694b967eb1SSatish Balay   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
13704b967eb1SSatish Balay   Mat a = aij->A, b = aij->B;
1371a008b906SSatish Balay   int ierr,s1,s2,s3;
1372a008b906SSatish Balay 
13734b967eb1SSatish Balay   ierr = MatGetLocalSize(mat,&s2,&s3); CHKERRQ(ierr);
13744b967eb1SSatish Balay   if (rr) {
13754b967eb1SSatish Balay     s3 = aij->n;
13764b967eb1SSatish Balay     VecGetLocalSize_Fast(rr,s1);
1377e3372554SBarry Smith     if (s1!=s3) SETERRQ(1,0,"right vector non-conforming local size");
13784b967eb1SSatish Balay     /* Overlap communication with computation. */
137943a90d84SBarry Smith     ierr = VecScatterBegin(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr);
1380a008b906SSatish Balay   }
13814b967eb1SSatish Balay   if (ll) {
13824b967eb1SSatish Balay     VecGetLocalSize_Fast(ll,s1);
1383e3372554SBarry Smith     if (s1!=s2) SETERRQ(1,0,"left vector non-conforming local size");
1384c351683dSSatish Balay     ierr = (*b->ops.diagonalscale)(b,ll,0); CHKERRQ(ierr);
13854b967eb1SSatish Balay   }
13864b967eb1SSatish Balay   /* scale  the diagonal block */
1387c351683dSSatish Balay   ierr = (*a->ops.diagonalscale)(a,ll,rr); CHKERRQ(ierr);
13884b967eb1SSatish Balay 
13894b967eb1SSatish Balay   if (rr) {
13904b967eb1SSatish Balay     /* Do a scatter end and then right scale the off-diagonal block */
139143a90d84SBarry Smith     ierr = VecScatterEnd(rr,aij->lvec,INSERT_VALUES,SCATTER_FORWARD,aij->Mvctx); CHKERRQ(ierr);
1392c351683dSSatish Balay     ierr = (*b->ops.diagonalscale)(b,0,aij->lvec); CHKERRQ(ierr);
13934b967eb1SSatish Balay   }
13944b967eb1SSatish Balay 
1395a008b906SSatish Balay   return 0;
1396a008b906SSatish Balay }
1397a008b906SSatish Balay 
1398a008b906SSatish Balay 
1399682d7d0cSBarry Smith extern int MatPrintHelp_SeqAIJ(Mat);
14005615d1e5SSatish Balay #undef __FUNC__
14015eea60f9SBarry Smith #define __FUNC__ "MatPrintHelp_MPIAIJ" /* ADIC Ignore */
14028f6be9afSLois Curfman McInnes int MatPrintHelp_MPIAIJ(Mat A)
1403682d7d0cSBarry Smith {
1404682d7d0cSBarry Smith   Mat_MPIAIJ *a   = (Mat_MPIAIJ*) A->data;
1405682d7d0cSBarry Smith 
1406682d7d0cSBarry Smith   if (!a->rank) return MatPrintHelp_SeqAIJ(a->A);
1407682d7d0cSBarry Smith   else return 0;
1408682d7d0cSBarry Smith }
1409682d7d0cSBarry Smith 
14105615d1e5SSatish Balay #undef __FUNC__
14115eea60f9SBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAIJ" /* ADIC Ignore */
14128f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIAIJ(Mat A,int *bs)
14135a838052SSatish Balay {
14145a838052SSatish Balay   *bs = 1;
14155a838052SSatish Balay   return 0;
14165a838052SSatish Balay }
14175615d1e5SSatish Balay #undef __FUNC__
14185eea60f9SBarry Smith #define __FUNC__ "MatSetUnfactored_MPIAIJ" /* ADIC Ignore */
14198f6be9afSLois Curfman McInnes int MatSetUnfactored_MPIAIJ(Mat A)
1420bb5a7306SBarry Smith {
1421bb5a7306SBarry Smith   Mat_MPIAIJ *a   = (Mat_MPIAIJ*) A->data;
1422bb5a7306SBarry Smith   int        ierr;
1423bb5a7306SBarry Smith   ierr = MatSetUnfactored(a->A); CHKERRQ(ierr);
1424bb5a7306SBarry Smith   return 0;
1425bb5a7306SBarry Smith }
1426bb5a7306SBarry Smith 
14278f6be9afSLois Curfman McInnes extern int MatConvertSameType_MPIAIJ(Mat,Mat *,int);
14282f86bd48SSatish Balay extern int MatIncreaseOverlap_MPIAIJ(Mat , int, IS *, int);
14290a198c4cSBarry Smith extern int MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring);
14300a198c4cSBarry Smith extern int MatGetSubMatrices_MPIAIJ (Mat ,int , IS *,IS *,MatGetSubMatrixCall,Mat **);
14318a729477SBarry Smith /* -------------------------------------------------------------------*/
14322ee70a88SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIAIJ,
143339e00950SLois Curfman McInnes        MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ,
143444a69424SLois Curfman McInnes        MatMult_MPIAIJ,MatMultAdd_MPIAIJ,
143544a69424SLois Curfman McInnes        MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ,
143636ce4990SBarry Smith        0,0,
143736ce4990SBarry Smith        0,0,
143836ce4990SBarry Smith        0,0,
143944a69424SLois Curfman McInnes        MatRelax_MPIAIJ,
1440b7c46309SBarry Smith        MatTranspose_MPIAIJ,
1441a66be287SLois Curfman McInnes        MatGetInfo_MPIAIJ,0,
1442a008b906SSatish Balay        MatGetDiagonal_MPIAIJ,MatDiagonalScale_MPIAIJ,MatNorm_MPIAIJ,
1443ee50ffe9SBarry Smith        MatAssemblyBegin_MPIAIJ,MatAssemblyEnd_MPIAIJ,
14441eb62cbbSBarry Smith        0,
1445299609e3SLois Curfman McInnes        MatSetOption_MPIAIJ,MatZeroEntries_MPIAIJ,MatZeroRows_MPIAIJ,
144636ce4990SBarry Smith        0,0,0,0,
1447d1710a03SLois Curfman McInnes        MatGetSize_MPIAIJ,MatGetLocalSize_MPIAIJ,MatGetOwnershipRange_MPIAIJ,
144836ce4990SBarry Smith        0,0,
144994a9d846SBarry Smith        0,0,MatConvertSameType_MPIAIJ,0,0,
1450b49de8d1SLois Curfman McInnes        0,0,0,
1451598137ffSSatish Balay        MatGetSubMatrices_MPIAIJ,MatIncreaseOverlap_MPIAIJ,MatGetValues_MPIAIJ,0,
1452052efed2SBarry Smith        MatPrintHelp_MPIAIJ,
14533b2fbd54SBarry Smith        MatScale_MPIAIJ,0,0,0,
14540a198c4cSBarry Smith        MatGetBlockSize_MPIAIJ,0,0,0,0,
1455bb5a7306SBarry Smith        MatFDColoringCreate_MPIAIJ,0,MatSetUnfactored_MPIAIJ};
145636ce4990SBarry Smith 
14578a729477SBarry Smith 
14585615d1e5SSatish Balay #undef __FUNC__
14595615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIAIJ"
14601987afe7SBarry Smith /*@C
1461ff756334SLois Curfman McInnes    MatCreateMPIAIJ - Creates a sparse parallel matrix in AIJ format
14623a511b96SLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
14633a511b96SLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameters
14643a511b96SLois Curfman McInnes    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
14653a511b96SLois Curfman McInnes    performance can be increased by more than a factor of 50.
14668a729477SBarry Smith 
14678a729477SBarry Smith    Input Parameters:
14681eb62cbbSBarry Smith .  comm - MPI communicator
14697d3e4905SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
147092e8d321SLois Curfman McInnes            This value should be the same as the local size used in creating the
147192e8d321SLois Curfman McInnes            y vector for the matrix-vector product y = Ax.
14721a3896d6SBarry Smith .  n - This value should be the same as the local size used in creating the
14731a3896d6SBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
14741a3896d6SBarry Smith        calculated if N is given) For square matrices n is almost always m.
14757d3e4905SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
147692e8d321SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
1477ab693e5aSLois Curfman McInnes .  d_nz - number of nonzeros per row in diagonal portion of local submatrix
1478ff756334SLois Curfman McInnes            (same for all local rows)
14792bd5e0b2SLois Curfman McInnes .  d_nzz - array containing the number of nonzeros in the various rows of the
148092e8d321SLois Curfman McInnes            diagonal portion of the local submatrix (possibly different for each row)
14812bd5e0b2SLois Curfman McInnes            or PETSC_NULL. You must leave room for the diagonal entry even if
14822bd5e0b2SLois Curfman McInnes            it is zero.
14832bd5e0b2SLois Curfman McInnes .  o_nz - number of nonzeros per row in the off-diagonal portion of local
1484ab693e5aSLois Curfman McInnes            submatrix (same for all local rows).
14852bd5e0b2SLois Curfman McInnes .  o_nzz - array containing the number of nonzeros in the various rows of the
14862bd5e0b2SLois Curfman McInnes            off-diagonal portion of the local submatrix (possibly different for
14872bd5e0b2SLois Curfman McInnes            each row) or PETSC_NULL.
14888a729477SBarry Smith 
1489ff756334SLois Curfman McInnes    Output Parameter:
149044cd7ae7SLois Curfman McInnes .  A - the matrix
14918a729477SBarry Smith 
1492ff756334SLois Curfman McInnes    Notes:
1493ff756334SLois Curfman McInnes    The AIJ format (also called the Yale sparse matrix format or
1494ff756334SLois Curfman McInnes    compressed row storage), is fully compatible with standard Fortran 77
14950002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
14960002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
1497ff756334SLois Curfman McInnes 
1498ff756334SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
1499ff756334SLois Curfman McInnes    (possibly both).
1500ff756334SLois Curfman McInnes 
15015511cfe3SLois Curfman McInnes    By default, this format uses inodes (identical nodes) when possible.
15025511cfe3SLois Curfman McInnes    We search for consecutive rows with the same nonzero structure, thereby
15035511cfe3SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
15045511cfe3SLois Curfman McInnes 
15055511cfe3SLois Curfman McInnes    Options Database Keys:
15065511cfe3SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
15076e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
15086e62573dSLois Curfman McInnes $        (max limit=5)
15096e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
15106e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
15116e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
15125511cfe3SLois Curfman McInnes 
1513e0245417SLois Curfman McInnes    Storage Information:
1514e0245417SLois Curfman McInnes    For a square global matrix we define each processor's diagonal portion
1515e0245417SLois Curfman McInnes    to be its local rows and the corresponding columns (a square submatrix);
1516e0245417SLois Curfman McInnes    each processor's off-diagonal portion encompasses the remainder of the
1517e0245417SLois Curfman McInnes    local matrix (a rectangular submatrix).
1518e0245417SLois Curfman McInnes 
1519e0245417SLois Curfman McInnes    The user can specify preallocated storage for the diagonal part of
15205ace5be8SLois Curfman McInnes    the local submatrix with either d_nz or d_nnz (not both).  Set
15215ace5be8SLois Curfman McInnes    d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
15225ace5be8SLois Curfman McInnes    memory allocation.  Likewise, specify preallocated storage for the
15235ace5be8SLois Curfman McInnes    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
1524ff756334SLois Curfman McInnes 
15255511cfe3SLois Curfman McInnes    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
15265511cfe3SLois Curfman McInnes    the figure below we depict these three local rows and all columns (0-11).
15272191d07cSBarry Smith 
1528b810aeb4SBarry Smith $          0 1 2 3 4 5 6 7 8 9 10 11
1529b810aeb4SBarry Smith $         -------------------
15305511cfe3SLois Curfman McInnes $  row 3  |  o o o d d d o o o o o o
15315511cfe3SLois Curfman McInnes $  row 4  |  o o o d d d o o o o o o
15325511cfe3SLois Curfman McInnes $  row 5  |  o o o d d d o o o o o o
15335511cfe3SLois Curfman McInnes $         -------------------
1534b810aeb4SBarry Smith $
1535b810aeb4SBarry Smith 
15365511cfe3SLois Curfman McInnes    Thus, any entries in the d locations are stored in the d (diagonal)
15375511cfe3SLois Curfman McInnes    submatrix, and any entries in the o locations are stored in the
15385511cfe3SLois Curfman McInnes    o (off-diagonal) submatrix.  Note that the d and the o submatrices are
15395511cfe3SLois Curfman McInnes    stored simply in the MATSEQAIJ format for compressed row storage.
15405511cfe3SLois Curfman McInnes 
15415511cfe3SLois Curfman McInnes    Now d_nz should indicate the number of nonzeros per row in the d matrix,
15425511cfe3SLois Curfman McInnes    and o_nz should indicate the number of nonzeros per row in the o matrix.
15435511cfe3SLois Curfman McInnes    In general, for PDE problems in which most nonzeros are near the diagonal,
15443d323bbdSBarry Smith    one expects d_nz >> o_nz. For large problems you MUST preallocate memory
154592e8d321SLois Curfman McInnes    or you will get TERRIBLE performance; see the users' manual chapter on
15466da5968aSLois Curfman McInnes    matrices.
15473a511b96SLois Curfman McInnes 
1548dbd7a890SLois Curfman McInnes .keywords: matrix, aij, compressed row, sparse, parallel
1549ff756334SLois Curfman McInnes 
1550fafbff53SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
15518a729477SBarry Smith @*/
15521eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N,
155344cd7ae7SLois Curfman McInnes                     int d_nz,int *d_nnz,int o_nz,int *o_nnz,Mat *A)
15548a729477SBarry Smith {
155544cd7ae7SLois Curfman McInnes   Mat          B;
155644cd7ae7SLois Curfman McInnes   Mat_MPIAIJ   *b;
155736ce4990SBarry Smith   int          ierr, i,sum[2],work[2],size;
1558416022c9SBarry Smith 
155944cd7ae7SLois Curfman McInnes   *A = 0;
156044cd7ae7SLois Curfman McInnes   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATMPIAIJ,comm);
156144cd7ae7SLois Curfman McInnes   PLogObjectCreate(B);
156244cd7ae7SLois Curfman McInnes   B->data       = (void *) (b = PetscNew(Mat_MPIAIJ)); CHKPTRQ(b);
156344cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_MPIAIJ));
156444cd7ae7SLois Curfman McInnes   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
156536ce4990SBarry Smith   MPI_Comm_size(comm,&size);
156636ce4990SBarry Smith   if (size == 1) {
156736ce4990SBarry Smith     B->ops.getrowij          = MatGetRowIJ_MPIAIJ;
156836ce4990SBarry Smith     B->ops.restorerowij      = MatRestoreRowIJ_MPIAIJ;
156936ce4990SBarry Smith     B->ops.lufactorsymbolic  = MatLUFactorSymbolic_MPIAIJ;
157036ce4990SBarry Smith     B->ops.lufactornumeric   = MatLUFactorNumeric_MPIAIJ;
157136ce4990SBarry Smith     B->ops.lufactor          = MatLUFactor_MPIAIJ;
157236ce4990SBarry Smith     B->ops.solve             = MatSolve_MPIAIJ;
157336ce4990SBarry Smith     B->ops.solveadd          = MatSolveAdd_MPIAIJ;
157436ce4990SBarry Smith     B->ops.solvetrans        = MatSolveTrans_MPIAIJ;
157536ce4990SBarry Smith     B->ops.solvetransadd     = MatSolveTransAdd_MPIAIJ;
157636ce4990SBarry Smith     B->ops.ilufactorsymbolic = MatILUFactorSymbolic_MPIAIJ;
157736ce4990SBarry Smith   }
157844cd7ae7SLois Curfman McInnes   B->destroy    = MatDestroy_MPIAIJ;
157944cd7ae7SLois Curfman McInnes   B->view       = MatView_MPIAIJ;
158044cd7ae7SLois Curfman McInnes   B->factor     = 0;
158144cd7ae7SLois Curfman McInnes   B->assembled  = PETSC_FALSE;
158290f02eecSBarry Smith   B->mapping    = 0;
1583d6dfbf8fSBarry Smith 
158447794344SBarry Smith   B->insertmode = NOT_SET_VALUES;
158544cd7ae7SLois Curfman McInnes   MPI_Comm_rank(comm,&b->rank);
158644cd7ae7SLois Curfman McInnes   MPI_Comm_size(comm,&b->size);
15871eb62cbbSBarry Smith 
1588b4fd4287SBarry Smith   if (m == PETSC_DECIDE && (d_nnz != PETSC_NULL || o_nnz != PETSC_NULL))
1589e3372554SBarry Smith     SETERRQ(1,0,"Cannot have PETSC_DECIDE rows but set d_nnz or o_nnz");
15901987afe7SBarry Smith 
1591dbd7a890SLois Curfman McInnes   if (M == PETSC_DECIDE || N == PETSC_DECIDE) {
15921eb62cbbSBarry Smith     work[0] = m; work[1] = n;
1593d65a2f8fSBarry Smith     MPI_Allreduce( work, sum,2,MPI_INT,MPI_SUM,comm );
1594dbd7a890SLois Curfman McInnes     if (M == PETSC_DECIDE) M = sum[0];
1595dbd7a890SLois Curfman McInnes     if (N == PETSC_DECIDE) N = sum[1];
15961eb62cbbSBarry Smith   }
159744cd7ae7SLois Curfman McInnes   if (m == PETSC_DECIDE) {m = M/b->size + ((M % b->size) > b->rank);}
159844cd7ae7SLois Curfman McInnes   if (n == PETSC_DECIDE) {n = N/b->size + ((N % b->size) > b->rank);}
159944cd7ae7SLois Curfman McInnes   b->m = m; B->m = m;
160044cd7ae7SLois Curfman McInnes   b->n = n; B->n = n;
160144cd7ae7SLois Curfman McInnes   b->N = N; B->N = N;
160244cd7ae7SLois Curfman McInnes   b->M = M; B->M = M;
16031eb62cbbSBarry Smith 
16041eb62cbbSBarry Smith   /* build local table of row and column ownerships */
160544cd7ae7SLois Curfman McInnes   b->rowners = (int *) PetscMalloc(2*(b->size+2)*sizeof(int)); CHKPTRQ(b->rowners);
160644cd7ae7SLois Curfman McInnes   PLogObjectMemory(B,2*(b->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ));
1607603f58a4SSatish Balay   b->cowners = b->rowners + b->size + 2;
160844cd7ae7SLois Curfman McInnes   MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);
160944cd7ae7SLois Curfman McInnes   b->rowners[0] = 0;
161044cd7ae7SLois Curfman McInnes   for ( i=2; i<=b->size; i++ ) {
161144cd7ae7SLois Curfman McInnes     b->rowners[i] += b->rowners[i-1];
16128a729477SBarry Smith   }
161344cd7ae7SLois Curfman McInnes   b->rstart = b->rowners[b->rank];
161444cd7ae7SLois Curfman McInnes   b->rend   = b->rowners[b->rank+1];
161544cd7ae7SLois Curfman McInnes   MPI_Allgather(&n,1,MPI_INT,b->cowners+1,1,MPI_INT,comm);
161644cd7ae7SLois Curfman McInnes   b->cowners[0] = 0;
161744cd7ae7SLois Curfman McInnes   for ( i=2; i<=b->size; i++ ) {
161844cd7ae7SLois Curfman McInnes     b->cowners[i] += b->cowners[i-1];
16198a729477SBarry Smith   }
162044cd7ae7SLois Curfman McInnes   b->cstart = b->cowners[b->rank];
162144cd7ae7SLois Curfman McInnes   b->cend   = b->cowners[b->rank+1];
16228a729477SBarry Smith 
16235ace5be8SLois Curfman McInnes   if (d_nz == PETSC_DEFAULT) d_nz = 5;
1624029af93fSBarry Smith   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,d_nz,d_nnz,&b->A); CHKERRQ(ierr);
162544cd7ae7SLois Curfman McInnes   PLogObjectParent(B,b->A);
16267b8455f0SLois Curfman McInnes   if (o_nz == PETSC_DEFAULT) o_nz = 0;
1627029af93fSBarry Smith   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m,N,o_nz,o_nnz,&b->B); CHKERRQ(ierr);
162844cd7ae7SLois Curfman McInnes   PLogObjectParent(B,b->B);
16298a729477SBarry Smith 
16301eb62cbbSBarry Smith   /* build cache for off array entries formed */
163144cd7ae7SLois Curfman McInnes   ierr = StashBuild_Private(&b->stash); CHKERRQ(ierr);
163290f02eecSBarry Smith   b->donotstash  = 0;
163344cd7ae7SLois Curfman McInnes   b->colmap      = 0;
163444cd7ae7SLois Curfman McInnes   b->garray      = 0;
163544cd7ae7SLois Curfman McInnes   b->roworiented = 1;
16368a729477SBarry Smith 
16371eb62cbbSBarry Smith   /* stuff used for matrix vector multiply */
163844cd7ae7SLois Curfman McInnes   b->lvec      = 0;
163944cd7ae7SLois Curfman McInnes   b->Mvctx     = 0;
16408a729477SBarry Smith 
16417a0afa10SBarry Smith   /* stuff for MatGetRow() */
164244cd7ae7SLois Curfman McInnes   b->rowindices   = 0;
164344cd7ae7SLois Curfman McInnes   b->rowvalues    = 0;
164444cd7ae7SLois Curfman McInnes   b->getrowactive = PETSC_FALSE;
16457a0afa10SBarry Smith 
164644cd7ae7SLois Curfman McInnes   *A = B;
1647d6dfbf8fSBarry Smith   return 0;
1648d6dfbf8fSBarry Smith }
1649c74985f6SBarry Smith 
16505615d1e5SSatish Balay #undef __FUNC__
16515615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIAIJ"
16528f6be9afSLois Curfman McInnes int MatConvertSameType_MPIAIJ(Mat matin,Mat *newmat,int cpvalues)
1653d6dfbf8fSBarry Smith {
1654d6dfbf8fSBarry Smith   Mat        mat;
1655416022c9SBarry Smith   Mat_MPIAIJ *a,*oldmat = (Mat_MPIAIJ *) matin->data;
1656a1b97e82SLois Curfman McInnes   int        ierr, len=0, flg;
1657d6dfbf8fSBarry Smith 
1658416022c9SBarry Smith   *newmat       = 0;
16590452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm);
1660a5a9c739SBarry Smith   PLogObjectCreate(mat);
16610452661fSBarry Smith   mat->data       = (void *) (a = PetscNew(Mat_MPIAIJ)); CHKPTRQ(a);
1662416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
166344a69424SLois Curfman McInnes   mat->destroy    = MatDestroy_MPIAIJ;
166444a69424SLois Curfman McInnes   mat->view       = MatView_MPIAIJ;
1665d6dfbf8fSBarry Smith   mat->factor     = matin->factor;
1666c456f294SBarry Smith   mat->assembled  = PETSC_TRUE;
1667d6dfbf8fSBarry Smith 
166844cd7ae7SLois Curfman McInnes   a->m = mat->m   = oldmat->m;
166944cd7ae7SLois Curfman McInnes   a->n = mat->n   = oldmat->n;
167044cd7ae7SLois Curfman McInnes   a->M = mat->M   = oldmat->M;
167144cd7ae7SLois Curfman McInnes   a->N = mat->N   = oldmat->N;
1672d6dfbf8fSBarry Smith 
1673416022c9SBarry Smith   a->rstart       = oldmat->rstart;
1674416022c9SBarry Smith   a->rend         = oldmat->rend;
1675416022c9SBarry Smith   a->cstart       = oldmat->cstart;
1676416022c9SBarry Smith   a->cend         = oldmat->cend;
167717699dbbSLois Curfman McInnes   a->size         = oldmat->size;
167817699dbbSLois Curfman McInnes   a->rank         = oldmat->rank;
167947794344SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1680bcd2baecSBarry Smith   a->rowvalues    = 0;
1681bcd2baecSBarry Smith   a->getrowactive = PETSC_FALSE;
1682d6dfbf8fSBarry Smith 
1683603f58a4SSatish Balay   a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners);
1684603f58a4SSatish Balay   PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIAIJ));
1685603f58a4SSatish Balay   a->cowners = a->rowners + a->size + 2;
1686603f58a4SSatish Balay   PetscMemcpy(a->rowners,oldmat->rowners,2*(a->size+2)*sizeof(int));
1687416022c9SBarry Smith   ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr);
16882ee70a88SLois Curfman McInnes   if (oldmat->colmap) {
16890452661fSBarry Smith     a->colmap = (int *) PetscMalloc((a->N)*sizeof(int));CHKPTRQ(a->colmap);
1690416022c9SBarry Smith     PLogObjectMemory(mat,(a->N)*sizeof(int));
1691416022c9SBarry Smith     PetscMemcpy(a->colmap,oldmat->colmap,(a->N)*sizeof(int));
1692416022c9SBarry Smith   } else a->colmap = 0;
16933f41c07dSBarry Smith   if (oldmat->garray) {
16943f41c07dSBarry Smith     len = ((Mat_SeqAIJ *) (oldmat->B->data))->n;
16953f41c07dSBarry Smith     a->garray = (int *) PetscMalloc((len+1)*sizeof(int)); CHKPTRQ(a->garray);
1696464493b3SBarry Smith     PLogObjectMemory(mat,len*sizeof(int));
16973f41c07dSBarry Smith     if (len) PetscMemcpy(a->garray,oldmat->garray,len*sizeof(int));
1698416022c9SBarry Smith   } else a->garray = 0;
1699d6dfbf8fSBarry Smith 
1700416022c9SBarry Smith   ierr =  VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr);
1701416022c9SBarry Smith   PLogObjectParent(mat,a->lvec);
1702a56f8943SBarry Smith   ierr =  VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr);
1703416022c9SBarry Smith   PLogObjectParent(mat,a->Mvctx);
1704416022c9SBarry Smith   ierr =  MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr);
1705416022c9SBarry Smith   PLogObjectParent(mat,a->A);
1706416022c9SBarry Smith   ierr =  MatConvert(oldmat->B,MATSAME,&a->B); CHKERRQ(ierr);
1707416022c9SBarry Smith   PLogObjectParent(mat,a->B);
17085dd7a6c7SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
170925cdf11fSBarry Smith   if (flg) {
1710682d7d0cSBarry Smith     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
1711682d7d0cSBarry Smith   }
17128a729477SBarry Smith   *newmat = mat;
17138a729477SBarry Smith   return 0;
17148a729477SBarry Smith }
1715416022c9SBarry Smith 
171677c4ece6SBarry Smith #include "sys.h"
1717416022c9SBarry Smith 
17185615d1e5SSatish Balay #undef __FUNC__
17195615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIAIJ"
172019bcc07fSBarry Smith int MatLoad_MPIAIJ(Viewer viewer,MatType type,Mat *newmat)
1721416022c9SBarry Smith {
1722d65a2f8fSBarry Smith   Mat          A;
1723416022c9SBarry Smith   int          i, nz, ierr, j,rstart, rend, fd;
1724d65a2f8fSBarry Smith   Scalar       *vals,*svals;
172519bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1726416022c9SBarry Smith   MPI_Status   status;
172717699dbbSLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
1728d65a2f8fSBarry Smith   int          *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols;
172919bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
1730416022c9SBarry Smith 
173117699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank);
173217699dbbSLois Curfman McInnes   if (!rank) {
173390ace30eSBarry Smith     ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
173477c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr);
1735e3372554SBarry Smith     if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object");
1736416022c9SBarry Smith   }
1737416022c9SBarry Smith 
1738416022c9SBarry Smith   MPI_Bcast(header+1,3,MPI_INT,0,comm);
1739416022c9SBarry Smith   M = header[1]; N = header[2];
1740416022c9SBarry Smith   /* determine ownership of all rows */
174117699dbbSLois Curfman McInnes   m = M/size + ((M % size) > rank);
17420452661fSBarry Smith   rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners);
1743416022c9SBarry Smith   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1744416022c9SBarry Smith   rowners[0] = 0;
174517699dbbSLois Curfman McInnes   for ( i=2; i<=size; i++ ) {
1746416022c9SBarry Smith     rowners[i] += rowners[i-1];
1747416022c9SBarry Smith   }
174817699dbbSLois Curfman McInnes   rstart = rowners[rank];
174917699dbbSLois Curfman McInnes   rend   = rowners[rank+1];
1750416022c9SBarry Smith 
1751416022c9SBarry Smith   /* distribute row lengths to all processors */
17520452661fSBarry Smith   ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens);
1753416022c9SBarry Smith   offlens = ourlens + (rend-rstart);
175417699dbbSLois Curfman McInnes   if (!rank) {
17550452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
175677c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
17570452661fSBarry Smith     sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts);
175817699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i];
1759d65a2f8fSBarry Smith     MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
17600452661fSBarry Smith     PetscFree(sndcounts);
1761416022c9SBarry Smith   }
1762416022c9SBarry Smith   else {
1763416022c9SBarry Smith     MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);
1764416022c9SBarry Smith   }
1765416022c9SBarry Smith 
176617699dbbSLois Curfman McInnes   if (!rank) {
1767416022c9SBarry Smith     /* calculate the number of nonzeros on each processor */
17680452661fSBarry Smith     procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz);
1769cddf8d76SBarry Smith     PetscMemzero(procsnz,size*sizeof(int));
177017699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) {
1771416022c9SBarry Smith       for ( j=rowners[i]; j< rowners[i+1]; j++ ) {
1772416022c9SBarry Smith         procsnz[i] += rowlengths[j];
1773416022c9SBarry Smith       }
1774416022c9SBarry Smith     }
17750452661fSBarry Smith     PetscFree(rowlengths);
1776416022c9SBarry Smith 
1777416022c9SBarry Smith     /* determine max buffer needed and allocate it */
1778416022c9SBarry Smith     maxnz = 0;
177917699dbbSLois Curfman McInnes     for ( i=0; i<size; i++ ) {
17800452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
1781416022c9SBarry Smith     }
17820452661fSBarry Smith     cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols);
1783416022c9SBarry Smith 
1784416022c9SBarry Smith     /* read in my part of the matrix column indices  */
1785416022c9SBarry Smith     nz = procsnz[0];
17860452661fSBarry Smith     mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
178777c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr);
1788d65a2f8fSBarry Smith 
1789d65a2f8fSBarry Smith     /* read in every one elses and ship off */
179017699dbbSLois Curfman McInnes     for ( i=1; i<size; i++ ) {
1791d65a2f8fSBarry Smith       nz = procsnz[i];
179277c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
1793d65a2f8fSBarry Smith       MPI_Send(cols,nz,MPI_INT,i,tag,comm);
1794d65a2f8fSBarry Smith     }
17950452661fSBarry Smith     PetscFree(cols);
1796416022c9SBarry Smith   }
1797416022c9SBarry Smith   else {
1798416022c9SBarry Smith     /* determine buffer space needed for message */
1799416022c9SBarry Smith     nz = 0;
1800416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
1801416022c9SBarry Smith       nz += ourlens[i];
1802416022c9SBarry Smith     }
18030452661fSBarry Smith     mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols);
1804416022c9SBarry Smith 
1805416022c9SBarry Smith     /* receive message of column indices*/
1806d65a2f8fSBarry Smith     MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);
1807416022c9SBarry Smith     MPI_Get_count(&status,MPI_INT,&maxnz);
1808e3372554SBarry Smith     if (maxnz != nz) SETERRQ(1,0,"something is wrong with file");
1809416022c9SBarry Smith   }
1810416022c9SBarry Smith 
1811416022c9SBarry Smith   /* loop over local rows, determining number of off diagonal entries */
1812cddf8d76SBarry Smith   PetscMemzero(offlens,m*sizeof(int));
1813416022c9SBarry Smith   jj = 0;
1814416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1815416022c9SBarry Smith     for ( j=0; j<ourlens[i]; j++ ) {
1816d65a2f8fSBarry Smith       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
1817416022c9SBarry Smith       jj++;
1818416022c9SBarry Smith     }
1819416022c9SBarry Smith   }
1820d65a2f8fSBarry Smith 
1821d65a2f8fSBarry Smith   /* create our matrix */
1822416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
1823416022c9SBarry Smith     ourlens[i] -= offlens[i];
1824416022c9SBarry Smith   }
1825d65a2f8fSBarry Smith   ierr = MatCreateMPIAIJ(comm,m,PETSC_DECIDE,M,N,0,ourlens,0,offlens,newmat);CHKERRQ(ierr);
1826d65a2f8fSBarry Smith   A = *newmat;
18276d4a8577SBarry Smith   MatSetOption(A,MAT_COLUMNS_SORTED);
1828d65a2f8fSBarry Smith   for ( i=0; i<m; i++ ) {
1829d65a2f8fSBarry Smith     ourlens[i] += offlens[i];
1830d65a2f8fSBarry Smith   }
1831416022c9SBarry Smith 
183217699dbbSLois Curfman McInnes   if (!rank) {
18330452661fSBarry Smith     vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals);
1834416022c9SBarry Smith 
1835416022c9SBarry Smith     /* read in my part of the matrix numerical values  */
1836416022c9SBarry Smith     nz = procsnz[0];
183777c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
1838d65a2f8fSBarry Smith 
1839d65a2f8fSBarry Smith     /* insert into matrix */
1840d65a2f8fSBarry Smith     jj      = rstart;
1841d65a2f8fSBarry Smith     smycols = mycols;
1842d65a2f8fSBarry Smith     svals   = vals;
1843d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1844d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1845d65a2f8fSBarry Smith       smycols += ourlens[i];
1846d65a2f8fSBarry Smith       svals   += ourlens[i];
1847d65a2f8fSBarry Smith       jj++;
1848416022c9SBarry Smith     }
1849416022c9SBarry Smith 
1850d65a2f8fSBarry Smith     /* read in other processors and ship out */
185117699dbbSLois Curfman McInnes     for ( i=1; i<size; i++ ) {
1852416022c9SBarry Smith       nz = procsnz[i];
185377c4ece6SBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
1854d65a2f8fSBarry Smith       MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
1855416022c9SBarry Smith     }
18560452661fSBarry Smith     PetscFree(procsnz);
1857416022c9SBarry Smith   }
1858d65a2f8fSBarry Smith   else {
1859d65a2f8fSBarry Smith     /* receive numeric values */
18600452661fSBarry Smith     vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
1861416022c9SBarry Smith 
1862d65a2f8fSBarry Smith     /* receive message of values*/
1863d65a2f8fSBarry Smith     MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
1864d65a2f8fSBarry Smith     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1865e3372554SBarry Smith     if (maxnz != nz) SETERRQ(1,0,"something is wrong with file");
1866d65a2f8fSBarry Smith 
1867d65a2f8fSBarry Smith     /* insert into matrix */
1868d65a2f8fSBarry Smith     jj      = rstart;
1869d65a2f8fSBarry Smith     smycols = mycols;
1870d65a2f8fSBarry Smith     svals   = vals;
1871d65a2f8fSBarry Smith     for ( i=0; i<m; i++ ) {
1872d65a2f8fSBarry Smith       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
1873d65a2f8fSBarry Smith       smycols += ourlens[i];
1874d65a2f8fSBarry Smith       svals   += ourlens[i];
1875d65a2f8fSBarry Smith       jj++;
1876d65a2f8fSBarry Smith     }
1877d65a2f8fSBarry Smith   }
18780452661fSBarry Smith   PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners);
1879d65a2f8fSBarry Smith 
18806d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
18816d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1882416022c9SBarry Smith   return 0;
1883416022c9SBarry Smith }
1884