xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 43ee02c34f45d7793d9e0dc9f8b4f1391ef178e7)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*43ee02c3SBarry Smith static char vcid[] = "$Id: aij.c,v 1.229 1997/07/09 20:53:48 balay Exp bsmith $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5d5d45c9bSBarry Smith /*
63369ce9aSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7d5d45c9bSBarry Smith   matrix storage format.
8d5d45c9bSBarry Smith */
93369ce9aSBarry Smith 
103369ce9aSBarry Smith #include "pinclude/pviewer.h"
113369ce9aSBarry Smith #include "sys.h"
1270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
13f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
14f5eb4b81SSatish Balay #include "src/inline/spops.h"
15f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1617ab2063SBarry Smith 
17a2ce50c7SBarry Smith /*
18a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
19a2ce50c7SBarry Smith */
205615d1e5SSatish Balay #undef __FUNC__
215615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ"
22a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
23a2ce50c7SBarry Smith {
24a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
25a2ce50c7SBarry Smith   int        ierr = 1;
26a2ce50c7SBarry Smith 
27e3372554SBarry Smith   SETERRQ(ierr,0,"Not implemented");
28a2ce50c7SBarry Smith }
29a2ce50c7SBarry Smith 
30bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3117ab2063SBarry Smith 
325615d1e5SSatish Balay #undef __FUNC__
335615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ"
348f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
356945ee14SBarry Smith                            PetscTruth *done)
3617ab2063SBarry Smith {
37416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
386945ee14SBarry Smith   int        ierr,i,ishift;
3917ab2063SBarry Smith 
4031625ec5SSatish Balay   *m     = A->m;
416945ee14SBarry Smith   if (!ia) return 0;
426945ee14SBarry Smith   ishift = a->indexshift;
436945ee14SBarry Smith   if (symmetric) {
4431625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
456945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
4631625ec5SSatish Balay     int nz = a->i[a->m];
473b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
4831625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
493b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
503b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
5131625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
526945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
5331625ec5SSatish Balay     int nz = a->i[a->m] + 1;
543b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
5531625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
563b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
573b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
5831625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
596945ee14SBarry Smith   } else {
606945ee14SBarry Smith     *ia = a->i; *ja = a->j;
61a2ce50c7SBarry Smith   }
62a2ce50c7SBarry Smith 
63a2744918SBarry Smith   return 0;
64a2744918SBarry Smith }
65a2744918SBarry Smith 
665615d1e5SSatish Balay #undef __FUNC__
675615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ"
688f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
696945ee14SBarry Smith                                PetscTruth *done)
706945ee14SBarry Smith {
716945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
723b2fbd54SBarry Smith   int        ishift = a->indexshift;
736945ee14SBarry Smith 
746945ee14SBarry Smith   if (!ia) return 0;
753b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
766945ee14SBarry Smith     PetscFree(*ia);
776945ee14SBarry Smith     PetscFree(*ja);
78bcd2baecSBarry Smith   }
7917ab2063SBarry Smith   return 0;
8017ab2063SBarry Smith }
8117ab2063SBarry Smith 
825615d1e5SSatish Balay #undef __FUNC__
835615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ"
8443a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
853b2fbd54SBarry Smith                            PetscTruth *done)
863b2fbd54SBarry Smith {
873b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
88a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
89a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
903b2fbd54SBarry Smith 
913b2fbd54SBarry Smith   *nn     = A->n;
923b2fbd54SBarry Smith   if (!ia) return 0;
933b2fbd54SBarry Smith   if (symmetric) {
94179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
953b2fbd54SBarry Smith   } else {
9661d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
973b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
983b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
99a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
1003b2fbd54SBarry Smith     jj = a->j;
1013b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
1023b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
1033b2fbd54SBarry Smith     }
1043b2fbd54SBarry Smith     cia[0] = oshift;
1053b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
1063b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
1073b2fbd54SBarry Smith     }
1083b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1093b2fbd54SBarry Smith     jj = a->j;
110a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
111a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
112a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1133b2fbd54SBarry Smith         col = *jj++ + ishift;
1143b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1153b2fbd54SBarry Smith       }
1163b2fbd54SBarry Smith     }
1173b2fbd54SBarry Smith     PetscFree(collengths);
1183b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1193b2fbd54SBarry Smith   }
1203b2fbd54SBarry Smith 
1213b2fbd54SBarry Smith   return 0;
1223b2fbd54SBarry Smith }
1233b2fbd54SBarry Smith 
1245615d1e5SSatish Balay #undef __FUNC__
1255615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ"
12643a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1273b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1283b2fbd54SBarry Smith {
1293b2fbd54SBarry Smith   if (!ia) return 0;
1303b2fbd54SBarry Smith 
1313b2fbd54SBarry Smith   PetscFree(*ia);
1323b2fbd54SBarry Smith   PetscFree(*ja);
1333b2fbd54SBarry Smith 
1343b2fbd54SBarry Smith   return 0;
1353b2fbd54SBarry Smith }
1363b2fbd54SBarry Smith 
137227d817aSBarry Smith #define CHUNKSIZE   15
13817ab2063SBarry Smith 
1395615d1e5SSatish Balay #undef __FUNC__
1405615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ"
14161d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
14217ab2063SBarry Smith {
143416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
144416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1454b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
146d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
147416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
14817ab2063SBarry Smith 
14917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
150416022c9SBarry Smith     row  = im[k];
1513b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
152e3372554SBarry Smith     if (row < 0) SETERRQ(1,0,"Negative row");
153e3372554SBarry Smith     if (row >= a->m) SETERRQ(1,0,"Row too large");
1543b2fbd54SBarry Smith #endif
15517ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
15617ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
157416022c9SBarry Smith     low = 0;
15817ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1593b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
160e3372554SBarry Smith       if (in[l] < 0) SETERRQ(1,0,"Negative column");
161e3372554SBarry Smith       if (in[l] >= a->n) SETERRQ(1,0,"Column too large");
1623b2fbd54SBarry Smith #endif
1634b0e389bSBarry Smith       col = in[l] - shift;
1644b0e389bSBarry Smith       if (roworiented) {
1654b0e389bSBarry Smith         value = *v++;
1664b0e389bSBarry Smith       }
1674b0e389bSBarry Smith       else {
1684b0e389bSBarry Smith         value = v[k + l*m];
1694b0e389bSBarry Smith       }
170416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
171416022c9SBarry Smith       while (high-low > 5) {
172416022c9SBarry Smith         t = (low+high)/2;
173416022c9SBarry Smith         if (rp[t] > col) high = t;
174416022c9SBarry Smith         else             low  = t;
17517ab2063SBarry Smith       }
176416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
17717ab2063SBarry Smith         if (rp[i] > col) break;
17817ab2063SBarry Smith         if (rp[i] == col) {
179416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
18017ab2063SBarry Smith           else                  ap[i] = value;
18117ab2063SBarry Smith           goto noinsert;
18217ab2063SBarry Smith         }
18317ab2063SBarry Smith       }
184c2653b3dSLois Curfman McInnes       if (nonew == 1) goto noinsert;
18511ebbc71SLois Curfman McInnes       else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix");
18617ab2063SBarry Smith       if (nrow >= rmax) {
18717ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
188416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
18917ab2063SBarry Smith         Scalar *new_a;
19017ab2063SBarry Smith 
19111ebbc71SLois Curfman McInnes         if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix");
19296854ed6SLois Curfman McInnes 
19317ab2063SBarry Smith         /* malloc new storage space */
194416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1950452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
19617ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
19717ab2063SBarry Smith         new_i   = new_j + new_nz;
19817ab2063SBarry Smith 
19917ab2063SBarry Smith         /* copy over old data into new slots */
20017ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
201416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
202416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
203416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
204416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
20517ab2063SBarry Smith                                                            len*sizeof(int));
206416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
207416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
20817ab2063SBarry Smith                                                            len*sizeof(Scalar));
20917ab2063SBarry Smith         /* free up old matrix storage */
2100452661fSBarry Smith         PetscFree(a->a);
2110452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
212416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
213416022c9SBarry Smith         a->singlemalloc = 1;
21417ab2063SBarry Smith 
21517ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
216416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
217416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
218416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
219b810aeb4SBarry Smith         a->reallocs++;
22017ab2063SBarry Smith       }
221416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
222416022c9SBarry Smith       /* shift up all the later entries in this row */
223416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
22417ab2063SBarry Smith         rp[ii+1] = rp[ii];
22517ab2063SBarry Smith         ap[ii+1] = ap[ii];
22617ab2063SBarry Smith       }
22717ab2063SBarry Smith       rp[i] = col;
22817ab2063SBarry Smith       ap[i] = value;
22917ab2063SBarry Smith       noinsert:;
230416022c9SBarry Smith       low = i + 1;
23117ab2063SBarry Smith     }
23217ab2063SBarry Smith     ailen[row] = nrow;
23317ab2063SBarry Smith   }
23417ab2063SBarry Smith   return 0;
23517ab2063SBarry Smith }
23617ab2063SBarry Smith 
2375615d1e5SSatish Balay #undef __FUNC__
2385615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ"
2398f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2407eb43aa7SLois Curfman McInnes {
2417eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
242b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2437eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2447eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2457eb43aa7SLois Curfman McInnes 
2467eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2477eb43aa7SLois Curfman McInnes     row  = im[k];
248e3372554SBarry Smith     if (row < 0) SETERRQ(1,0,"Negative row");
249e3372554SBarry Smith     if (row >= a->m) SETERRQ(1,0,"Row too large");
2507eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2517eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2527eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
253e3372554SBarry Smith       if (in[l] < 0) SETERRQ(1,0,"Negative column");
254e3372554SBarry Smith       if (in[l] >= a->n) SETERRQ(1,0,"Column too large");
2557eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2567eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2577eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2587eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2597eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2607eb43aa7SLois Curfman McInnes         else             low  = t;
2617eb43aa7SLois Curfman McInnes       }
2627eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2637eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2647eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
265b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2667eb43aa7SLois Curfman McInnes           goto finished;
2677eb43aa7SLois Curfman McInnes         }
2687eb43aa7SLois Curfman McInnes       }
269b49de8d1SLois Curfman McInnes       *v++ = zero;
2707eb43aa7SLois Curfman McInnes       finished:;
2717eb43aa7SLois Curfman McInnes     }
2727eb43aa7SLois Curfman McInnes   }
2737eb43aa7SLois Curfman McInnes   return 0;
2747eb43aa7SLois Curfman McInnes }
2757eb43aa7SLois Curfman McInnes 
27617ab2063SBarry Smith 
2775615d1e5SSatish Balay #undef __FUNC__
2785615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary"
2798f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
28017ab2063SBarry Smith {
281416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
282416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
28317ab2063SBarry Smith 
28490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2850452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
286416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
287416022c9SBarry Smith   col_lens[1] = a->m;
288416022c9SBarry Smith   col_lens[2] = a->n;
289416022c9SBarry Smith   col_lens[3] = a->nz;
290416022c9SBarry Smith 
291416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
292416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
293416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
29417ab2063SBarry Smith   }
29577c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
2960452661fSBarry Smith   PetscFree(col_lens);
297416022c9SBarry Smith 
298416022c9SBarry Smith   /* store column indices (zero start index) */
299416022c9SBarry Smith   if (a->indexshift) {
300416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
30117ab2063SBarry Smith   }
30277c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr);
303416022c9SBarry Smith   if (a->indexshift) {
304416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
30517ab2063SBarry Smith   }
306416022c9SBarry Smith 
307416022c9SBarry Smith   /* store nonzero values */
30877c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
30917ab2063SBarry Smith   return 0;
31017ab2063SBarry Smith }
311416022c9SBarry Smith 
3125615d1e5SSatish Balay #undef __FUNC__
3135615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII"
3148f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
315416022c9SBarry Smith {
316416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
317496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
31817ab2063SBarry Smith   FILE        *fd;
31917ab2063SBarry Smith   char        *outputname;
32017ab2063SBarry Smith 
32190ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
322416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
32390ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
324a93ec695SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
32595e01e2fSLois Curfman McInnes     return 0;
32695e01e2fSLois Curfman McInnes   }
327a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
328496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
329496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
330496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
33195e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
33295e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
33317ab2063SBarry Smith   }
334a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
335416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3364e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
3374e220ebcSLois Curfman McInnes     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz);
33817ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
33917ab2063SBarry Smith 
34017ab2063SBarry Smith     for (i=0; i<m; i++) {
341416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
34217ab2063SBarry Smith #if defined(PETSC_COMPLEX)
3436945ee14SBarry Smith         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]),
344416022c9SBarry Smith                    imag(a->a[j]));
34517ab2063SBarry Smith #else
3467a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
34717ab2063SBarry Smith #endif
34817ab2063SBarry Smith       }
34917ab2063SBarry Smith     }
35017ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
35117ab2063SBarry Smith   }
352a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
35344cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
35444cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
35544cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
35644cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
357766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0 && real(a->a[j]) != 0.0)
35844cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
359766eeae4SLois Curfman McInnes         else if (imag(a->a[j]) < 0.0 && real(a->a[j]) != 0.0)
360766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
36144cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
36244cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
36344cd7ae7SLois Curfman McInnes #else
36444cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
36544cd7ae7SLois Curfman McInnes #endif
36644cd7ae7SLois Curfman McInnes       }
36744cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
36844cd7ae7SLois Curfman McInnes     }
36944cd7ae7SLois Curfman McInnes   }
370496be53dSLois Curfman McInnes   else if (format == VIEWER_FORMAT_ASCII_SYMMODU) {
371496be53dSLois Curfman McInnes     int nzd=0, fshift=1, *sptr;
3722e44a96cSLois Curfman McInnes     sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr);
373496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
374496be53dSLois Curfman McInnes       sptr[i] = nzd+1;
375496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
376496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
377496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX)
378496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) nzd++;
379496be53dSLois Curfman McInnes #else
380496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) nzd++;
381496be53dSLois Curfman McInnes #endif
382496be53dSLois Curfman McInnes         }
383496be53dSLois Curfman McInnes       }
384496be53dSLois Curfman McInnes     }
3852e44a96cSLois Curfman McInnes     sptr[m] = nzd+1;
386496be53dSLois Curfman McInnes     fprintf(fd," %d %d\n\n",m,nzd);
3872e44a96cSLois Curfman McInnes     for ( i=0; i<m+1; i+=6 ) {
3882e44a96cSLois Curfman McInnes       if (i+4<m) fprintf(fd," %d %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);
3892e44a96cSLois Curfman McInnes       else if (i+3<m) fprintf(fd," %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);
3902e44a96cSLois Curfman McInnes       else if (i+2<m) fprintf(fd," %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);
3912e44a96cSLois Curfman McInnes       else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]);
3922e44a96cSLois Curfman McInnes       else if (i<m)   fprintf(fd," %d %d\n",sptr[i],sptr[i+1]);
3937272d637SLois Curfman McInnes       else            fprintf(fd," %d\n",sptr[i]);
394496be53dSLois Curfman McInnes     }
395496be53dSLois Curfman McInnes     fprintf(fd,"\n");
396496be53dSLois Curfman McInnes     PetscFree(sptr);
397496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
398496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
399496be53dSLois Curfman McInnes         if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift);
400496be53dSLois Curfman McInnes       }
401496be53dSLois Curfman McInnes       fprintf(fd,"\n");
402496be53dSLois Curfman McInnes     }
403496be53dSLois Curfman McInnes     fprintf(fd,"\n");
404496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
405496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
406496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
407496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX)
408496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0)
409496be53dSLois Curfman McInnes             fprintf(fd," %18.16e %18.16e ",real(a->a[j]),imag(a->a[j]));
410496be53dSLois Curfman McInnes #else
411496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]);
412496be53dSLois Curfman McInnes #endif
413496be53dSLois Curfman McInnes         }
414496be53dSLois Curfman McInnes       }
415496be53dSLois Curfman McInnes       fprintf(fd,"\n");
416496be53dSLois Curfman McInnes     }
417496be53dSLois Curfman McInnes   }
41817ab2063SBarry Smith   else {
41917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
42017ab2063SBarry Smith       fprintf(fd,"row %d:",i);
421416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
42217ab2063SBarry Smith #if defined(PETSC_COMPLEX)
423766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0) {
424416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
425766eeae4SLois Curfman McInnes         } else if (imag(a->a[j]) < 0.0) {
426766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
42717ab2063SBarry Smith         }
42817ab2063SBarry Smith         else {
429416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
43017ab2063SBarry Smith         }
43117ab2063SBarry Smith #else
432416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
43317ab2063SBarry Smith #endif
43417ab2063SBarry Smith       }
43517ab2063SBarry Smith       fprintf(fd,"\n");
43617ab2063SBarry Smith     }
43717ab2063SBarry Smith   }
43817ab2063SBarry Smith   fflush(fd);
439416022c9SBarry Smith   return 0;
440416022c9SBarry Smith }
441416022c9SBarry Smith 
4425615d1e5SSatish Balay #undef __FUNC__
4435615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Draw"
4448f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
445416022c9SBarry Smith {
446416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
447cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
4480513a670SBarry Smith   int         format;
44994a9d846SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r,maxv = 0.0;
450bcd2baecSBarry Smith   Draw        draw;
451cddf8d76SBarry Smith   DrawButton  button;
45219bcc07fSBarry Smith   PetscTruth  isnull;
453cddf8d76SBarry Smith 
4540513a670SBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
4550513a670SBarry Smith   ierr = DrawClear(draw); CHKERRQ(ierr);
4560513a670SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
45719bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
45819bcc07fSBarry Smith 
459416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
460416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
461416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
462416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
4630513a670SBarry Smith 
4640513a670SBarry Smith   if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4650513a670SBarry Smith     /* Blue for negative, Cyan for zero and  Red for positive */
466cddf8d76SBarry Smith     color = DRAW_BLUE;
467416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
468cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
469416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
470cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
471cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
472cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
473cddf8d76SBarry Smith #else
474cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
475cddf8d76SBarry Smith #endif
476cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
477cddf8d76SBarry Smith       }
478cddf8d76SBarry Smith     }
479cddf8d76SBarry Smith     color = DRAW_CYAN;
480cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
481cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
482cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
483cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
484cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
485cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
486cddf8d76SBarry Smith       }
487cddf8d76SBarry Smith     }
488cddf8d76SBarry Smith     color = DRAW_RED;
489cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
490cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
491cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
492cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
493cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
494cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
495cddf8d76SBarry Smith #else
496cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
497cddf8d76SBarry Smith #endif
498cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
499416022c9SBarry Smith       }
500416022c9SBarry Smith     }
5010513a670SBarry Smith   } else {
5020513a670SBarry Smith     /* use contour shading to indicate magnitude of values */
5030513a670SBarry Smith     /* first determine max of all nonzero values */
5040513a670SBarry Smith     int    nz = a->nz,count;
5050513a670SBarry Smith     Draw   popup;
5060513a670SBarry Smith 
5070513a670SBarry Smith     for ( i=0; i<nz; i++ ) {
5080513a670SBarry Smith       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
5090513a670SBarry Smith     }
5100513a670SBarry Smith     ierr = DrawCreatePopUp(draw,&popup); CHKERRQ(ierr);
5110513a670SBarry Smith     ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr);
5120513a670SBarry Smith     count = 0;
5130513a670SBarry Smith     for ( i=0; i<m; i++ ) {
5140513a670SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
5150513a670SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5160513a670SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5170513a670SBarry Smith         color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5180513a670SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5190513a670SBarry Smith         count++;
5200513a670SBarry Smith       }
5210513a670SBarry Smith     }
5220513a670SBarry Smith   }
523416022c9SBarry Smith   DrawFlush(draw);
524cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
525cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
526cddf8d76SBarry Smith 
527cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
5286945ee14SBarry Smith   ierr = DrawCheckResizedWindow(draw);
529cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
530cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
531cddf8d76SBarry Smith     DrawClear(draw);
532cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
533cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
534cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
535cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
536cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
537cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
538cddf8d76SBarry Smith     w *= scale; h *= scale;
539cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5400513a670SBarry Smith     if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
5410513a670SBarry Smith       /* Blue for negative, Cyan for zero and  Red for positive */
542cddf8d76SBarry Smith       color = DRAW_BLUE;
543cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
544cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
545cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
546cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
547cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
548cddf8d76SBarry Smith           if (real(a->a[j]) >=  0.) continue;
549cddf8d76SBarry Smith #else
550cddf8d76SBarry Smith           if (a->a[j] >=  0.) continue;
551cddf8d76SBarry Smith #endif
552cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
553cddf8d76SBarry Smith         }
554cddf8d76SBarry Smith       }
555cddf8d76SBarry Smith       color = DRAW_CYAN;
556cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
557cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
558cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
559cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
560cddf8d76SBarry Smith           if (a->a[j] !=  0.) continue;
561cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
562cddf8d76SBarry Smith         }
563cddf8d76SBarry Smith       }
564cddf8d76SBarry Smith       color = DRAW_RED;
565cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
566cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
567cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
568cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
569cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
570cddf8d76SBarry Smith           if (real(a->a[j]) <=  0.) continue;
571cddf8d76SBarry Smith #else
572cddf8d76SBarry Smith           if (a->a[j] <=  0.) continue;
573cddf8d76SBarry Smith #endif
574cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
575cddf8d76SBarry Smith         }
576cddf8d76SBarry Smith       }
5770513a670SBarry Smith     } else {
5780513a670SBarry Smith       /* use contour shading to indicate magnitude of values */
5790513a670SBarry Smith       int count = 0;
5800513a670SBarry Smith       for ( i=0; i<m; i++ ) {
5810513a670SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
5820513a670SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5830513a670SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
5840513a670SBarry Smith           color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5850513a670SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5860513a670SBarry Smith           count++;
5870513a670SBarry Smith         }
5880513a670SBarry Smith       }
5890513a670SBarry Smith     }
5900513a670SBarry Smith 
5916945ee14SBarry Smith     ierr = DrawCheckResizedWindow(draw);
592cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
593cddf8d76SBarry Smith   }
594416022c9SBarry Smith   return 0;
595416022c9SBarry Smith }
596416022c9SBarry Smith 
5975615d1e5SSatish Balay #undef __FUNC__
598c22c1629SBarry Smith #define __FUNC__ "MatView_SeqAIJ" /* ADIC Ignore */
5998f6be9afSLois Curfman McInnes int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
600416022c9SBarry Smith {
601416022c9SBarry Smith   Mat         A = (Mat) obj;
602416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
603bcd2baecSBarry Smith   ViewerType  vtype;
604bcd2baecSBarry Smith   int         ierr;
605416022c9SBarry Smith 
606bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
607bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
608416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
609416022c9SBarry Smith   }
610bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
611416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
612416022c9SBarry Smith   }
613bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
614416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
615416022c9SBarry Smith   }
616bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
617bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
61817ab2063SBarry Smith   }
61917ab2063SBarry Smith   return 0;
62017ab2063SBarry Smith }
62119bcc07fSBarry Smith 
622c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
6235615d1e5SSatish Balay #undef __FUNC__
6245615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ"
6258f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
62617ab2063SBarry Smith {
627416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
62841c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
629*43ee02c3SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0;
630416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
63117ab2063SBarry Smith 
6326d4a8577SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
63317ab2063SBarry Smith 
634*43ee02c3SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
63517ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
636416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
63717ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
63894a9d846SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
63917ab2063SBarry Smith     if (fshift) {
640416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
64117ab2063SBarry Smith       N = ailen[i];
64217ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
64317ab2063SBarry Smith         ip[j-fshift] = ip[j];
64417ab2063SBarry Smith         ap[j-fshift] = ap[j];
64517ab2063SBarry Smith       }
64617ab2063SBarry Smith     }
64717ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
64817ab2063SBarry Smith   }
64917ab2063SBarry Smith   if (m) {
65017ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
65117ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
65217ab2063SBarry Smith   }
65317ab2063SBarry Smith   /* reset ilen and imax for each row */
65417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
65517ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
65617ab2063SBarry Smith   }
657416022c9SBarry Smith   a->nz = ai[m] + shift;
65817ab2063SBarry Smith 
65917ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
660416022c9SBarry Smith   if (fshift && a->diag) {
6610452661fSBarry Smith     PetscFree(a->diag);
662416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
663416022c9SBarry Smith     a->diag = 0;
66417ab2063SBarry Smith   }
6654e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
6664e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
6674e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
668b810aeb4SBarry Smith            a->reallocs);
66994a9d846SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax);
670dd5f02e7SSatish Balay   a->reallocs          = 0;
6714e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
6724e220ebcSLois Curfman McInnes 
67376dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
67441c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
67517ab2063SBarry Smith   return 0;
67617ab2063SBarry Smith }
67717ab2063SBarry Smith 
6785615d1e5SSatish Balay #undef __FUNC__
6795615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ"
6808f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A)
68117ab2063SBarry Smith {
682416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
683cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
68417ab2063SBarry Smith   return 0;
68517ab2063SBarry Smith }
686416022c9SBarry Smith 
6875615d1e5SSatish Balay #undef __FUNC__
6885615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ"
68917ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
69017ab2063SBarry Smith {
691416022c9SBarry Smith   Mat        A  = (Mat) obj;
692416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
693d5d45c9bSBarry Smith 
69417ab2063SBarry Smith #if defined(PETSC_LOG)
695416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
69617ab2063SBarry Smith #endif
6970452661fSBarry Smith   PetscFree(a->a);
6980452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
6990452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
7000452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
7010452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
7020452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
70376dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
7040452661fSBarry Smith   PetscFree(a);
705eed86810SBarry Smith 
706f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
707f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
70817ab2063SBarry Smith   return 0;
70917ab2063SBarry Smith }
71017ab2063SBarry Smith 
7115615d1e5SSatish Balay #undef __FUNC__
7125615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
7138f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A)
71417ab2063SBarry Smith {
71517ab2063SBarry Smith   return 0;
71617ab2063SBarry Smith }
71717ab2063SBarry Smith 
7185615d1e5SSatish Balay #undef __FUNC__
7195615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
7208f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op)
72117ab2063SBarry Smith {
722416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7236d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
7246d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
7256d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
726219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
7276d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
728c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
72996854ed6SLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
7306d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
7316d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
732219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7336d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7346d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
73590f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
7362b362799SSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES)
73794a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
7386d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
739e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");}
7406d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
7416d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
7426d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
7436d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
7446d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
745e2f28af5SBarry Smith   else
746e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
74717ab2063SBarry Smith   return 0;
74817ab2063SBarry Smith }
74917ab2063SBarry Smith 
7505615d1e5SSatish Balay #undef __FUNC__
7515615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
7528f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
75317ab2063SBarry Smith {
754416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
755416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
75617ab2063SBarry Smith   Scalar     *x, zero = 0.0;
75717ab2063SBarry Smith 
75817ab2063SBarry Smith   VecSet(&zero,v);
75990f02eecSBarry Smith   VecGetArray_Fast(v,x); VecGetLocalSize(v,&n);
760e3372554SBarry Smith   if (n != a->m) SETERRQ(1,0,"Nonconforming matrix and vector");
761416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
762416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
763416022c9SBarry Smith       if (a->j[j]+shift == i) {
764416022c9SBarry Smith         x[i] = a->a[j];
76517ab2063SBarry Smith         break;
76617ab2063SBarry Smith       }
76717ab2063SBarry Smith     }
76817ab2063SBarry Smith   }
76917ab2063SBarry Smith   return 0;
77017ab2063SBarry Smith }
77117ab2063SBarry Smith 
77217ab2063SBarry Smith /* -------------------------------------------------------*/
77317ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
77417ab2063SBarry Smith /* -------------------------------------------------------*/
7755615d1e5SSatish Balay #undef __FUNC__
7765615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
77744cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
77817ab2063SBarry Smith {
779416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
78017ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
781416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
78217ab2063SBarry Smith 
78390f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
784cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
78517ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
78617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
787416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
788416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
789416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
79017ab2063SBarry Smith     alpha = x[i];
79117ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
79217ab2063SBarry Smith   }
793416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
79417ab2063SBarry Smith   return 0;
79517ab2063SBarry Smith }
796d5d45c9bSBarry Smith 
7975615d1e5SSatish Balay #undef __FUNC__
7985615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
79944cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
80017ab2063SBarry Smith {
801416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
80217ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
803416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
80417ab2063SBarry Smith 
80590f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
80617ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
80717ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
80817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
809416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
810416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
811416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
81217ab2063SBarry Smith     alpha = x[i];
81317ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
81417ab2063SBarry Smith   }
81590f02eecSBarry Smith   PLogFlops(2*a->nz);
81617ab2063SBarry Smith   return 0;
81717ab2063SBarry Smith }
81817ab2063SBarry Smith 
8195615d1e5SSatish Balay #undef __FUNC__
8205615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
82144cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
82217ab2063SBarry Smith {
823416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
82417ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
8259ea0dfa2SSatish Balay   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j;
82617ab2063SBarry Smith 
82790f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
82817ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
829416022c9SBarry Smith   idx  = a->j;
830416022c9SBarry Smith   v    = a->a;
831416022c9SBarry Smith   ii   = a->i;
83217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8339ea0dfa2SSatish Balay     jrow = ii[i];
8349ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
83517ab2063SBarry Smith     sum  = 0.0;
8369ea0dfa2SSatish Balay     /* while (n--) sum += *v++ * x[*idx++]; */
8379ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8389ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8399ea0dfa2SSatish Balay      }
84017ab2063SBarry Smith     y[i] = sum;
84117ab2063SBarry Smith   }
842416022c9SBarry Smith   PLogFlops(2*a->nz - m);
84317ab2063SBarry Smith   return 0;
84417ab2063SBarry Smith }
84517ab2063SBarry Smith 
8465615d1e5SSatish Balay #undef __FUNC__
8475615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
84844cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
84917ab2063SBarry Smith {
850416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
85117ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
8529ea0dfa2SSatish Balay   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j;
8539ea0dfa2SSatish Balay 
85417ab2063SBarry Smith 
85590f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); VecGetArray_Fast(zz,z);
85617ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
857cddf8d76SBarry Smith   idx  = a->j;
858cddf8d76SBarry Smith   v    = a->a;
859cddf8d76SBarry Smith   ii   = a->i;
86017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8619ea0dfa2SSatish Balay     jrow = ii[i];
8629ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
86317ab2063SBarry Smith     sum  = y[i];
8649ea0dfa2SSatish Balay     /* while (n--) sum += *v++ * x[*idx++]; */
8659ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8669ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8679ea0dfa2SSatish Balay      }
86817ab2063SBarry Smith     z[i] = sum;
86917ab2063SBarry Smith   }
870416022c9SBarry Smith   PLogFlops(2*a->nz);
87117ab2063SBarry Smith   return 0;
87217ab2063SBarry Smith }
87317ab2063SBarry Smith 
87417ab2063SBarry Smith /*
87517ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
87617ab2063SBarry Smith */
87717ab2063SBarry Smith 
8785615d1e5SSatish Balay #undef __FUNC__
8795615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
880416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
88117ab2063SBarry Smith {
882416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
883416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
88417ab2063SBarry Smith 
8850452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
886416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
887416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
888416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
889416022c9SBarry Smith       if (a->j[j]+shift == i) {
89017ab2063SBarry Smith         diag[i] = j - shift;
89117ab2063SBarry Smith         break;
89217ab2063SBarry Smith       }
89317ab2063SBarry Smith     }
89417ab2063SBarry Smith   }
895416022c9SBarry Smith   a->diag = diag;
89617ab2063SBarry Smith   return 0;
89717ab2063SBarry Smith }
89817ab2063SBarry Smith 
8995615d1e5SSatish Balay #undef __FUNC__
9005615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
90144cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
90217ab2063SBarry Smith                            double fshift,int its,Vec xx)
90317ab2063SBarry Smith {
904416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
905416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
906d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
90717ab2063SBarry Smith 
90890f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(bb,b);
909416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
910416022c9SBarry Smith   diag = a->diag;
911416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
91217ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
91317ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
91417ab2063SBarry Smith     bs = b + shift;
91517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
916416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
917416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
918416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
919416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
92017ab2063SBarry Smith         sum  = b[i]*d/omega;
92117ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
92217ab2063SBarry Smith         x[i] = sum;
92317ab2063SBarry Smith     }
92417ab2063SBarry Smith     return 0;
92517ab2063SBarry Smith   }
92617ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
927e3372554SBarry Smith     SETERRQ(1,0,"SOR_APPLY_LOWER is not done");
92817ab2063SBarry Smith   }
929416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
93017ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
93117ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
93217ab2063SBarry Smith 
93317ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
93417ab2063SBarry Smith 
93517ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
93617ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
93717ab2063SBarry Smith     is the relaxation factor.
93817ab2063SBarry Smith     */
9390452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
94017ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
94117ab2063SBarry Smith 
94217ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
94317ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
944416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
945416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
946416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
947416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
94817ab2063SBarry Smith       sum  = b[i];
94917ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
95017ab2063SBarry Smith       x[i] = omega*(sum/d);
95117ab2063SBarry Smith     }
95217ab2063SBarry Smith 
95317ab2063SBarry Smith     /*  t = b - (2*E - D)x */
954416022c9SBarry Smith     v = a->a;
95517ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
95617ab2063SBarry Smith 
95717ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
958416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
959416022c9SBarry Smith     diag = a->diag;
96017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
961416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
962416022c9SBarry Smith       n    = diag[i] - a->i[i];
963416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
964416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
96517ab2063SBarry Smith       sum  = t[i];
96617ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
96717ab2063SBarry Smith       t[i] = omega*(sum/d);
96817ab2063SBarry Smith     }
96917ab2063SBarry Smith 
97017ab2063SBarry Smith     /*  x = x + t */
97117ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
9720452661fSBarry Smith     PetscFree(t);
97317ab2063SBarry Smith     return 0;
97417ab2063SBarry Smith   }
97517ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
97617ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
97717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
978416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
979416022c9SBarry Smith         n    = diag[i] - a->i[i];
980416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
981416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
98217ab2063SBarry Smith         sum  = b[i];
98317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
98417ab2063SBarry Smith         x[i] = omega*(sum/d);
98517ab2063SBarry Smith       }
98617ab2063SBarry Smith       xb = x;
98717ab2063SBarry Smith     }
98817ab2063SBarry Smith     else xb = b;
98917ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
99017ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
99117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
992416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
99317ab2063SBarry Smith       }
99417ab2063SBarry Smith     }
99517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
99617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
997416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
998416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
999416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1000416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
100117ab2063SBarry Smith         sum  = xb[i];
100217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
100317ab2063SBarry Smith         x[i] = omega*(sum/d);
100417ab2063SBarry Smith       }
100517ab2063SBarry Smith     }
100617ab2063SBarry Smith     its--;
100717ab2063SBarry Smith   }
100817ab2063SBarry Smith   while (its--) {
100917ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
101017ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1011416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1012416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1013416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1014416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
101517ab2063SBarry Smith         sum  = b[i];
101617ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10177e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
101817ab2063SBarry Smith       }
101917ab2063SBarry Smith     }
102017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
102117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1022416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1023416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1024416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1025416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
102617ab2063SBarry Smith         sum  = b[i];
102717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10287e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
102917ab2063SBarry Smith       }
103017ab2063SBarry Smith     }
103117ab2063SBarry Smith   }
103217ab2063SBarry Smith   return 0;
103317ab2063SBarry Smith }
103417ab2063SBarry Smith 
10355615d1e5SSatish Balay #undef __FUNC__
10365615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10378f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
103817ab2063SBarry Smith {
1039416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
10404e220ebcSLois Curfman McInnes 
10414e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
10424e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
10434e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
10444e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
10454e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10464e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
10474e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
10484e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
10494e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
10504e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
10514e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
10524e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
10534e220ebcSLois Curfman McInnes   info->memory         = A->mem;
10544e220ebcSLois Curfman McInnes   if (A->factor) {
10554e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10564e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
10574e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
10584e220ebcSLois Curfman McInnes   } else {
10594e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
10604e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
10614e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
10624e220ebcSLois Curfman McInnes   }
106317ab2063SBarry Smith   return 0;
106417ab2063SBarry Smith }
106517ab2063SBarry Smith 
106617ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
106717ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
106817ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
106917ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
107017ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
107117ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
107217ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
107317ab2063SBarry Smith 
10745615d1e5SSatish Balay #undef __FUNC__
10755615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
10768f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
107717ab2063SBarry Smith {
1078416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1079416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
108017ab2063SBarry Smith 
108177c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
108217ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
108317ab2063SBarry Smith   if (diag) {
108417ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1085e3372554SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range");
1086416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1087416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1088416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1089416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
109017ab2063SBarry Smith       }
109117ab2063SBarry Smith       else {
109217ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
109317ab2063SBarry Smith         CHKERRQ(ierr);
109417ab2063SBarry Smith       }
109517ab2063SBarry Smith     }
109617ab2063SBarry Smith   }
109717ab2063SBarry Smith   else {
109817ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1099e3372554SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range");
1100416022c9SBarry Smith       a->ilen[rows[i]] = 0;
110117ab2063SBarry Smith     }
110217ab2063SBarry Smith   }
110317ab2063SBarry Smith   ISRestoreIndices(is,&rows);
110443a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
110517ab2063SBarry Smith   return 0;
110617ab2063SBarry Smith }
110717ab2063SBarry Smith 
11085615d1e5SSatish Balay #undef __FUNC__
11095615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11108f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
111117ab2063SBarry Smith {
1112416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1113416022c9SBarry Smith   *m = a->m; *n = a->n;
111417ab2063SBarry Smith   return 0;
111517ab2063SBarry Smith }
111617ab2063SBarry Smith 
11175615d1e5SSatish Balay #undef __FUNC__
11185615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11198f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
112017ab2063SBarry Smith {
1121416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1122416022c9SBarry Smith   *m = 0; *n = a->m;
112317ab2063SBarry Smith   return 0;
112417ab2063SBarry Smith }
1125026e39d0SSatish Balay 
11265615d1e5SSatish Balay #undef __FUNC__
11275615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11284e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
112917ab2063SBarry Smith {
1130416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1131c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
113217ab2063SBarry Smith 
1133e3372554SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,0,"Row out of range");
113417ab2063SBarry Smith 
1135416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1136416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
113717ab2063SBarry Smith   if (idx) {
1138416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
11394e093b46SBarry Smith     if (*nz && shift) {
11400452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
114117ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
11424e093b46SBarry Smith     } else if (*nz) {
11434e093b46SBarry Smith       *idx = itmp;
114417ab2063SBarry Smith     }
114517ab2063SBarry Smith     else *idx = 0;
114617ab2063SBarry Smith   }
114717ab2063SBarry Smith   return 0;
114817ab2063SBarry Smith }
114917ab2063SBarry Smith 
11505615d1e5SSatish Balay #undef __FUNC__
11515615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
11524e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
115317ab2063SBarry Smith {
11544e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11554e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
115617ab2063SBarry Smith   return 0;
115717ab2063SBarry Smith }
115817ab2063SBarry Smith 
11595615d1e5SSatish Balay #undef __FUNC__
11605615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
11618f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
116217ab2063SBarry Smith {
1163416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1164416022c9SBarry Smith   Scalar     *v = a->a;
116517ab2063SBarry Smith   double     sum = 0.0;
1166416022c9SBarry Smith   int        i, j,shift = a->indexshift;
116717ab2063SBarry Smith 
116817ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1169416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
117017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
117117ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
117217ab2063SBarry Smith #else
117317ab2063SBarry Smith       sum += (*v)*(*v); v++;
117417ab2063SBarry Smith #endif
117517ab2063SBarry Smith     }
117617ab2063SBarry Smith     *norm = sqrt(sum);
117717ab2063SBarry Smith   }
117817ab2063SBarry Smith   else if (type == NORM_1) {
117917ab2063SBarry Smith     double *tmp;
1180416022c9SBarry Smith     int    *jj = a->j;
118166963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1182cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
118317ab2063SBarry Smith     *norm = 0.0;
1184416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1185a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
118617ab2063SBarry Smith     }
1187416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
118817ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
118917ab2063SBarry Smith     }
11900452661fSBarry Smith     PetscFree(tmp);
119117ab2063SBarry Smith   }
119217ab2063SBarry Smith   else if (type == NORM_INFINITY) {
119317ab2063SBarry Smith     *norm = 0.0;
1194416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1195416022c9SBarry Smith       v = a->a + a->i[j] + shift;
119617ab2063SBarry Smith       sum = 0.0;
1197416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1198cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
119917ab2063SBarry Smith       }
120017ab2063SBarry Smith       if (sum > *norm) *norm = sum;
120117ab2063SBarry Smith     }
120217ab2063SBarry Smith   }
120317ab2063SBarry Smith   else {
1204e3372554SBarry Smith     SETERRQ(1,0,"No support for two norm yet");
120517ab2063SBarry Smith   }
120617ab2063SBarry Smith   return 0;
120717ab2063SBarry Smith }
120817ab2063SBarry Smith 
12095615d1e5SSatish Balay #undef __FUNC__
12105615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12118f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
121217ab2063SBarry Smith {
1213416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1214416022c9SBarry Smith   Mat        C;
1215416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1216416022c9SBarry Smith   int        shift = a->indexshift;
1217d5d45c9bSBarry Smith   Scalar     *array = a->a;
121817ab2063SBarry Smith 
12193638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
1220e3372554SBarry Smith     SETERRQ(1,0,"Square matrix only for in-place");
12210452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1222cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
122317ab2063SBarry Smith   if (shift) {
122417ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
122517ab2063SBarry Smith   }
122617ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1227416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12280452661fSBarry Smith   PetscFree(col);
122917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
123017ab2063SBarry Smith     len = ai[i+1]-ai[i];
1231416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
123217ab2063SBarry Smith     array += len; aj += len;
123317ab2063SBarry Smith   }
123417ab2063SBarry Smith   if (shift) {
123517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
123617ab2063SBarry Smith   }
123717ab2063SBarry Smith 
12386d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12396d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
124017ab2063SBarry Smith 
12413638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1242416022c9SBarry Smith     *B = C;
124317ab2063SBarry Smith   } else {
1244416022c9SBarry Smith     /* This isn't really an in-place transpose */
12450452661fSBarry Smith     PetscFree(a->a);
12460452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
12470452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
12480452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
12490452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
12500452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
12511073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
12520452661fSBarry Smith     PetscFree(a);
1253f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
12540452661fSBarry Smith     PetscHeaderDestroy(C);
125517ab2063SBarry Smith   }
125617ab2063SBarry Smith   return 0;
125717ab2063SBarry Smith }
125817ab2063SBarry Smith 
12595615d1e5SSatish Balay #undef __FUNC__
12605615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
12618f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
126217ab2063SBarry Smith {
1263416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
126417ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1265d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
126617ab2063SBarry Smith 
126717ab2063SBarry Smith   if (ll) {
12683ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
12693ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
12709b1297e1SSatish Balay     VecGetLocalSize_Fast(ll,m);
1271e3372554SBarry Smith     if (m != a->m) SETERRQ(1,0,"Left scaling vector wrong length");
127290f02eecSBarry Smith     VecGetArray_Fast(ll,l);
1273416022c9SBarry Smith     v = a->a;
127417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
127517ab2063SBarry Smith       x = l[i];
1276416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
127717ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
127817ab2063SBarry Smith     }
127944cd7ae7SLois Curfman McInnes     PLogFlops(nz);
128017ab2063SBarry Smith   }
128117ab2063SBarry Smith   if (rr) {
12829b1297e1SSatish Balay     VecGetLocalSize_Fast(rr,n);
1283e3372554SBarry Smith     if (n != a->n) SETERRQ(1,0,"Right scaling vector wrong length");
128490f02eecSBarry Smith     VecGetArray_Fast(rr,r);
1285416022c9SBarry Smith     v = a->a; jj = a->j;
128617ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
128717ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
128817ab2063SBarry Smith     }
128944cd7ae7SLois Curfman McInnes     PLogFlops(nz);
129017ab2063SBarry Smith   }
129117ab2063SBarry Smith   return 0;
129217ab2063SBarry Smith }
129317ab2063SBarry Smith 
12945615d1e5SSatish Balay #undef __FUNC__
12955615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
12968f6be9afSLois Curfman McInnes int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
129717ab2063SBarry Smith {
1298db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
129902834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
130099141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1301a2744918SBarry Smith   register int sum,lensi;
130299141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
130399141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
130499141d43SSatish Balay   Scalar       *a_new,*mat_a;
1305416022c9SBarry Smith   Mat          C;
130617ab2063SBarry Smith 
1307b48a1e75SSatish Balay   ierr = ISSorted(isrow,(PetscTruth*)&i);
1308e3372554SBarry Smith   if (!i) SETERRQ(1,0,"ISrow is not sorted");
130999141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
1310e3372554SBarry Smith   if (!i) SETERRQ(1,0,"IScol is not sorted");
131199141d43SSatish Balay 
131217ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
131317ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
131417ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
131517ab2063SBarry Smith 
13167264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
131702834360SBarry Smith     /* special case of contiguous rows */
131857faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
131902834360SBarry Smith     starts = lens + ncols;
132002834360SBarry Smith     /* loop over new rows determining lens and starting points */
132102834360SBarry Smith     for (i=0; i<nrows; i++) {
1322a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1323a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
132402834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1325d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
132602834360SBarry Smith           starts[i] = k;
132702834360SBarry Smith           break;
132802834360SBarry Smith 	}
132902834360SBarry Smith       }
1330a2744918SBarry Smith       sum = 0;
133102834360SBarry Smith       while (k < kend) {
1332d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1333a2744918SBarry Smith         sum++;
133402834360SBarry Smith       }
1335a2744918SBarry Smith       lens[i] = sum;
133602834360SBarry Smith     }
133702834360SBarry Smith     /* create submatrix */
1338cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
133908480c60SBarry Smith       int n_cols,n_rows;
134008480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1341e3372554SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,0,"");
1342d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
134308480c60SBarry Smith       C = *B;
134408480c60SBarry Smith     }
134508480c60SBarry Smith     else {
134602834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
134708480c60SBarry Smith     }
1348db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1349db02288aSLois Curfman McInnes 
135002834360SBarry Smith     /* loop over rows inserting into submatrix */
1351db02288aSLois Curfman McInnes     a_new    = c->a;
1352db02288aSLois Curfman McInnes     j_new    = c->j;
1353db02288aSLois Curfman McInnes     i_new    = c->i;
1354db02288aSLois Curfman McInnes     i_new[0] = -shift;
135502834360SBarry Smith     for (i=0; i<nrows; i++) {
1356a2744918SBarry Smith       ii    = starts[i];
1357a2744918SBarry Smith       lensi = lens[i];
1358a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1359a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
136002834360SBarry Smith       }
1361a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1362a2744918SBarry Smith       a_new      += lensi;
1363a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1364a2744918SBarry Smith       c->ilen[i]  = lensi;
136502834360SBarry Smith     }
13660452661fSBarry Smith     PetscFree(lens);
136702834360SBarry Smith   }
136802834360SBarry Smith   else {
136902834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
13700452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
137102834360SBarry Smith     ssmap = smap + shift;
137299141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1373cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
137417ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
137502834360SBarry Smith     /* determine lens of each row */
137602834360SBarry Smith     for (i=0; i<nrows; i++) {
1377d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
137802834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
137902834360SBarry Smith       lens[i] = 0;
138002834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1381d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
138202834360SBarry Smith           lens[i]++;
138302834360SBarry Smith         }
138402834360SBarry Smith       }
138502834360SBarry Smith     }
138617ab2063SBarry Smith     /* Create and fill new matrix */
1387a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
138899141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
138999141d43SSatish Balay 
1390e3372554SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(1,0,"");
139199141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1392e3372554SBarry Smith         SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros");
139399141d43SSatish Balay       }
139499141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
139508480c60SBarry Smith       C = *B;
139608480c60SBarry Smith     }
139708480c60SBarry Smith     else {
139802834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
139908480c60SBarry Smith     }
140099141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
140117ab2063SBarry Smith     for (i=0; i<nrows; i++) {
140299141d43SSatish Balay       row    = irow[i];
140317ab2063SBarry Smith       nznew  = 0;
140499141d43SSatish Balay       kstart = ai[row]+shift;
140599141d43SSatish Balay       kend   = kstart + a->ilen[row];
140699141d43SSatish Balay       mat_i  = c->i[i]+shift;
140799141d43SSatish Balay       mat_j  = c->j + mat_i;
140899141d43SSatish Balay       mat_a  = c->a + mat_i;
140999141d43SSatish Balay       mat_ilen = c->ilen + i;
141017ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
141199141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
141299141d43SSatish Balay           *mat_j++ = tcol - (!shift);
141399141d43SSatish Balay           *mat_a++ = a->a[k];
141499141d43SSatish Balay           (*mat_ilen)++;
141599141d43SSatish Balay 
141617ab2063SBarry Smith         }
141717ab2063SBarry Smith       }
141817ab2063SBarry Smith     }
141902834360SBarry Smith     /* Free work space */
142002834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
142199141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
142202834360SBarry Smith   }
14236d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
14246d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
142517ab2063SBarry Smith 
142617ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1427416022c9SBarry Smith   *B = C;
142817ab2063SBarry Smith   return 0;
142917ab2063SBarry Smith }
143017ab2063SBarry Smith 
1431a871dcd8SBarry Smith /*
143263b91edcSBarry Smith      note: This can only work for identity for row and col. It would
143363b91edcSBarry Smith    be good to check this and otherwise generate an error.
1434a871dcd8SBarry Smith */
14355615d1e5SSatish Balay #undef __FUNC__
14365615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
14378f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1438a871dcd8SBarry Smith {
143963b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
144008480c60SBarry Smith   int        ierr;
144163b91edcSBarry Smith   Mat        outA;
144263b91edcSBarry Smith 
1443e3372554SBarry Smith   if (fill != 0) SETERRQ(1,0,"Only fill=0 supported");
1444a871dcd8SBarry Smith 
144563b91edcSBarry Smith   outA          = inA;
144663b91edcSBarry Smith   inA->factor   = FACTOR_LU;
144763b91edcSBarry Smith   a->row        = row;
144863b91edcSBarry Smith   a->col        = col;
144963b91edcSBarry Smith 
145094a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
14510452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
145294a9d846SBarry Smith   }
145363b91edcSBarry Smith 
145408480c60SBarry Smith   if (!a->diag) {
145508480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
145663b91edcSBarry Smith   }
145763b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1458a871dcd8SBarry Smith   return 0;
1459a871dcd8SBarry Smith }
1460a871dcd8SBarry Smith 
1461f0b747eeSBarry Smith #include "pinclude/plapack.h"
14625615d1e5SSatish Balay #undef __FUNC__
14635615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
14648f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1465f0b747eeSBarry Smith {
1466f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1467f0b747eeSBarry Smith   int        one = 1;
1468f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1469f0b747eeSBarry Smith   PLogFlops(a->nz);
1470f0b747eeSBarry Smith   return 0;
1471f0b747eeSBarry Smith }
1472f0b747eeSBarry Smith 
14735615d1e5SSatish Balay #undef __FUNC__
14745615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
14758f6be9afSLois Curfman McInnes int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1476cddf8d76SBarry Smith                                     Mat **B)
1477cddf8d76SBarry Smith {
1478cddf8d76SBarry Smith   int ierr,i;
1479cddf8d76SBarry Smith 
1480cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
14810452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1482cddf8d76SBarry Smith   }
1483cddf8d76SBarry Smith 
1484cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1485905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1486cddf8d76SBarry Smith   }
1487cddf8d76SBarry Smith   return 0;
1488cddf8d76SBarry Smith }
1489cddf8d76SBarry Smith 
14905615d1e5SSatish Balay #undef __FUNC__
14915615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
14928f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
14935a838052SSatish Balay {
14945a838052SSatish Balay   *bs = 1;
14955a838052SSatish Balay   return 0;
14965a838052SSatish Balay }
14975a838052SSatish Balay 
14985615d1e5SSatish Balay #undef __FUNC__
14995615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15008f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15014dcbc457SBarry Smith {
1502e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
150306763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15048a047759SSatish Balay   int        start, end, *ai, *aj;
150506763907SSatish Balay   char       *table;
15068a047759SSatish Balay   shift = a->indexshift;
1507e4d965acSSatish Balay   m     = a->m;
1508e4d965acSSatish Balay   ai    = a->i;
15098a047759SSatish Balay   aj    = a->j+shift;
15108a047759SSatish Balay 
1511e3372554SBarry Smith   if (ov < 0)  SETERRQ(1,0,"illegal overlap value used");
151206763907SSatish Balay 
151306763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
151406763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
151506763907SSatish Balay 
1516e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1517b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1518e4d965acSSatish Balay     isz  = 0;
151906763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1520e4d965acSSatish Balay 
1521e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
15224dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
152377c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1524e4d965acSSatish Balay 
1525dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1526e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
152706763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
15284dcbc457SBarry Smith     }
152906763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
153006763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1531e4d965acSSatish Balay 
153204a348a9SBarry Smith     k = 0;
153304a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
153404a348a9SBarry Smith       n = isz;
153506763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1536e4d965acSSatish Balay         row   = nidx[k];
1537e4d965acSSatish Balay         start = ai[row];
1538e4d965acSSatish Balay         end   = ai[row+1];
153904a348a9SBarry Smith         for ( l = start; l<end ; l++){
15408a047759SSatish Balay           val = aj[l] + shift;
154106763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1542e4d965acSSatish Balay         }
1543e4d965acSSatish Balay       }
1544e4d965acSSatish Balay     }
1545029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1546e4d965acSSatish Balay   }
154704a348a9SBarry Smith   PetscFree(table);
154806763907SSatish Balay   PetscFree(nidx);
1549e4d965acSSatish Balay   return 0;
15504dcbc457SBarry Smith }
155117ab2063SBarry Smith 
15520513a670SBarry Smith /* -------------------------------------------------------------- */
15530513a670SBarry Smith #undef __FUNC__
15540513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
15550513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
15560513a670SBarry Smith {
15570513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
15580513a670SBarry Smith   Scalar     *vwork;
15590513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
15600513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
156156cd22aeSBarry Smith   IS         icolp,irowp;
15620513a670SBarry Smith 
156356cd22aeSBarry Smith   ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr);
156456cd22aeSBarry Smith   ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr);
156556cd22aeSBarry Smith   ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr);
156656cd22aeSBarry Smith   ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr);
15670513a670SBarry Smith 
15680513a670SBarry Smith   /* determine lengths of permuted rows */
15690513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
15700513a670SBarry Smith   for (i=0; i<m; i++ ) {
15710513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
15720513a670SBarry Smith   }
15730513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
15740513a670SBarry Smith   PetscFree(lens);
15750513a670SBarry Smith 
15760513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
15770513a670SBarry Smith   for (i=0; i<m; i++) {
15780513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15790513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
15800513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
15810513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15820513a670SBarry Smith   }
15830513a670SBarry Smith   PetscFree(cnew);
15840513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
15850513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
158656cd22aeSBarry Smith   ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr);
158756cd22aeSBarry Smith   ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr);
158856cd22aeSBarry Smith   ierr = ISDestroy(irowp); CHKERRQ(ierr);
158956cd22aeSBarry Smith   ierr = ISDestroy(icolp); CHKERRQ(ierr);
15900513a670SBarry Smith   return 0;
15910513a670SBarry Smith }
15920513a670SBarry Smith 
15935615d1e5SSatish Balay #undef __FUNC__
15945615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1595682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1596682d7d0cSBarry Smith {
1597682d7d0cSBarry Smith   static int called = 0;
1598682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1599682d7d0cSBarry Smith 
1600682d7d0cSBarry Smith   if (called) return 0; else called = 1;
160177c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
16020f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
16030f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
16040f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_no_inode: Do not use inodes\n");
16050f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1606682d7d0cSBarry Smith #if defined(HAVE_ESSL)
16070f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1608682d7d0cSBarry Smith #endif
1609682d7d0cSBarry Smith   return 0;
1610682d7d0cSBarry Smith }
16118f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1612a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1613a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1614a93ec695SBarry Smith 
1615682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
161617ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
161717ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1618416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1619416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
162017ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
162117ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
162217ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
162317ab2063SBarry Smith        MatRelax_SeqAIJ,
162417ab2063SBarry Smith        MatTranspose_SeqAIJ,
16257264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1626f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
162717ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
162817ab2063SBarry Smith        MatCompress_SeqAIJ,
162917ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
163017ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
163117ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
163217ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
163394a9d846SBarry Smith        0,0,
16343d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1635cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
16367eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1637682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1638f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
16395a838052SSatish Balay        MatScale_SeqAIJ,0,0,
16406945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
16416945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
16423b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
16433b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
16443b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1645a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1646a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
16470513a670SBarry Smith        MatColoringPatch_SeqAIJ,
16480513a670SBarry Smith        0,
16490513a670SBarry Smith        MatPermute_SeqAIJ};
165017ab2063SBarry Smith 
165117ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
165217ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
165317ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
165417ab2063SBarry Smith 
16555615d1e5SSatish Balay #undef __FUNC__
16565615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
165717ab2063SBarry Smith /*@C
1658682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
16590d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
16606e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
16612bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
16622bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
166317ab2063SBarry Smith 
166417ab2063SBarry Smith    Input Parameters:
1665029af93fSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
166617ab2063SBarry Smith .  m - number of rows
166717ab2063SBarry Smith .  n - number of columns
166817ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
16692bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
16702bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
167117ab2063SBarry Smith 
167217ab2063SBarry Smith    Output Parameter:
1673416022c9SBarry Smith .  A - the matrix
167417ab2063SBarry Smith 
167517ab2063SBarry Smith    Notes:
167617ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
167717ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
16780002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
167944cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
168017ab2063SBarry Smith 
168117ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1682a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
16833d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
16846da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
168517ab2063SBarry Smith 
1686682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1687682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1688682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
16896c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
16906c7ebb05SLois Curfman McInnes 
16916c7ebb05SLois Curfman McInnes    Options Database Keys:
16926c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
16936e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
16946e62573dSLois Curfman McInnes $        (max limit=5)
16956e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
16966e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
16976e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
169817ab2063SBarry Smith 
169917ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
170017ab2063SBarry Smith @*/
1701416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
170217ab2063SBarry Smith {
1703416022c9SBarry Smith   Mat        B;
1704416022c9SBarry Smith   Mat_SeqAIJ *b;
17056945ee14SBarry Smith   int        i, len, ierr, flg,size;
17066945ee14SBarry Smith 
17076945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1708e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
1709d5d45c9bSBarry Smith 
1710416022c9SBarry Smith   *A                  = 0;
1711f09e8eb9SSatish Balay   PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1712416022c9SBarry Smith   PLogObjectCreate(B);
17130452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
171444cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1715416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1716416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1717416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1718416022c9SBarry Smith   B->factor           = 0;
1719416022c9SBarry Smith   B->lupivotthreshold = 1.0;
172090f02eecSBarry Smith   B->mapping          = 0;
17217a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
172269957df2SSatish Balay                           &flg); CHKERRQ(ierr);
17237a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
17247a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
17257a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1726416022c9SBarry Smith   b->row              = 0;
1727416022c9SBarry Smith   b->col              = 0;
1728416022c9SBarry Smith   b->indexshift       = 0;
1729b810aeb4SBarry Smith   b->reallocs         = 0;
173069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
173169957df2SSatish Balay   if (flg) b->indexshift = -1;
173217ab2063SBarry Smith 
173344cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
173444cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
17350452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1736b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
17377b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
17387b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1739416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
174017ab2063SBarry Smith     nz = nz*m;
174117ab2063SBarry Smith   }
174217ab2063SBarry Smith   else {
174317ab2063SBarry Smith     nz = 0;
1744416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
174517ab2063SBarry Smith   }
174617ab2063SBarry Smith 
174717ab2063SBarry Smith   /* allocate the matrix space */
1748416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
17490452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1750416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1751cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1752416022c9SBarry Smith   b->i  = b->j + nz;
1753416022c9SBarry Smith   b->singlemalloc = 1;
175417ab2063SBarry Smith 
1755416022c9SBarry Smith   b->i[0] = -b->indexshift;
175617ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1757416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
175817ab2063SBarry Smith   }
175917ab2063SBarry Smith 
1760416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
17610452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1762f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1763416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
176417ab2063SBarry Smith 
1765416022c9SBarry Smith   b->nz               = 0;
1766416022c9SBarry Smith   b->maxnz            = nz;
1767416022c9SBarry Smith   b->sorted           = 0;
1768416022c9SBarry Smith   b->roworiented      = 1;
1769416022c9SBarry Smith   b->nonew            = 0;
1770416022c9SBarry Smith   b->diag             = 0;
1771416022c9SBarry Smith   b->solve_work       = 0;
177271bd300dSLois Curfman McInnes   b->spptr            = 0;
1773754ec7b1SSatish Balay   b->inode.node_count = 0;
1774754ec7b1SSatish Balay   b->inode.size       = 0;
17756c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
17766c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
17774e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
177817ab2063SBarry Smith 
1779416022c9SBarry Smith   *A = B;
17804e220ebcSLois Curfman McInnes 
17814b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
17824b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
178369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
178469957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
17854b14c69eSBarry Smith #endif
178669957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
178769957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
178869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
178969957df2SSatish Balay   if (flg) {
1790e3372554SBarry Smith     if (!b->indexshift) SETERRQ(1,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1791416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
179217ab2063SBarry Smith   }
179369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
179469957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
179517ab2063SBarry Smith   return 0;
179617ab2063SBarry Smith }
179717ab2063SBarry Smith 
17985615d1e5SSatish Balay #undef __FUNC__
17995615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
18003d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
180117ab2063SBarry Smith {
1802416022c9SBarry Smith   Mat        C;
1803416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
180408480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
180517ab2063SBarry Smith 
18064043dd9cSLois Curfman McInnes   *B = 0;
1807f09e8eb9SSatish Balay   PetscHeaderCreate(C,_p_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1808416022c9SBarry Smith   PLogObjectCreate(C);
18090452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
181041c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1811416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1812416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1813416022c9SBarry Smith   C->factor     = A->factor;
1814416022c9SBarry Smith   c->row        = 0;
1815416022c9SBarry Smith   c->col        = 0;
1816416022c9SBarry Smith   c->indexshift = shift;
1817c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
181817ab2063SBarry Smith 
181944cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
182044cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
182144cd7ae7SLois Curfman McInnes   C->M          = a->m;
182244cd7ae7SLois Curfman McInnes   C->N          = a->n;
182317ab2063SBarry Smith 
18240452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
18250452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
182617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1827416022c9SBarry Smith     c->imax[i] = a->imax[i];
1828416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
182917ab2063SBarry Smith   }
183017ab2063SBarry Smith 
183117ab2063SBarry Smith   /* allocate the matrix space */
1832416022c9SBarry Smith   c->singlemalloc = 1;
1833416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
18340452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1835416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1836416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1837416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
183817ab2063SBarry Smith   if (m > 0) {
1839416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
184008480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1841416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
184217ab2063SBarry Smith     }
184308480c60SBarry Smith   }
184417ab2063SBarry Smith 
1845f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1846416022c9SBarry Smith   c->sorted      = a->sorted;
1847416022c9SBarry Smith   c->roworiented = a->roworiented;
1848416022c9SBarry Smith   c->nonew       = a->nonew;
18497a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
185017ab2063SBarry Smith 
1851416022c9SBarry Smith   if (a->diag) {
18520452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1853416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
185417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1855416022c9SBarry Smith       c->diag[i] = a->diag[i];
185617ab2063SBarry Smith     }
185717ab2063SBarry Smith   }
1858416022c9SBarry Smith   else c->diag          = 0;
18596c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
18606c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1861754ec7b1SSatish Balay   if (a->inode.size){
1862daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1863754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1864daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1865754ec7b1SSatish Balay   } else {
1866754ec7b1SSatish Balay     c->inode.size       = 0;
1867754ec7b1SSatish Balay     c->inode.node_count = 0;
1868754ec7b1SSatish Balay   }
1869416022c9SBarry Smith   c->nz                 = a->nz;
1870416022c9SBarry Smith   c->maxnz              = a->maxnz;
1871416022c9SBarry Smith   c->solve_work         = 0;
187276dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1873754ec7b1SSatish Balay 
1874416022c9SBarry Smith   *B = C;
187517ab2063SBarry Smith   return 0;
187617ab2063SBarry Smith }
187717ab2063SBarry Smith 
18785615d1e5SSatish Balay #undef __FUNC__
18795615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
188019bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
188117ab2063SBarry Smith {
1882416022c9SBarry Smith   Mat_SeqAIJ   *a;
1883416022c9SBarry Smith   Mat          B;
188417699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1885bcd2baecSBarry Smith   MPI_Comm     comm;
188617ab2063SBarry Smith 
188719bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
188817699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
1889e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"view must have one processor");
189090ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
189177c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1892e3372554SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object in file");
189317ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
189417ab2063SBarry Smith 
189517ab2063SBarry Smith   /* read in row lengths */
18960452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
189777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
189817ab2063SBarry Smith 
189917ab2063SBarry Smith   /* create our matrix */
1900416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1901416022c9SBarry Smith   B = *A;
1902416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1903416022c9SBarry Smith   shift = a->indexshift;
190417ab2063SBarry Smith 
190517ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
190677c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
190717ab2063SBarry Smith   if (shift) {
190817ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1909416022c9SBarry Smith       a->j[i] += 1;
191017ab2063SBarry Smith     }
191117ab2063SBarry Smith   }
191217ab2063SBarry Smith 
191317ab2063SBarry Smith   /* read in nonzero values */
191477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
191517ab2063SBarry Smith 
191617ab2063SBarry Smith   /* set matrix "i" values */
1917416022c9SBarry Smith   a->i[0] = -shift;
191817ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1919416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1920416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
192117ab2063SBarry Smith   }
19220452661fSBarry Smith   PetscFree(rowlengths);
192317ab2063SBarry Smith 
19246d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
19256d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
192617ab2063SBarry Smith   return 0;
192717ab2063SBarry Smith }
192817ab2063SBarry Smith 
19295615d1e5SSatish Balay #undef __FUNC__
19305615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
19318f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
19327264ac53SSatish Balay {
19337264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
19347264ac53SSatish Balay 
1935e3372554SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,0,"Matrices must be same type");
19367264ac53SSatish Balay 
19377264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
19387264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1939bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
194077c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1941bcd2baecSBarry Smith   }
19427264ac53SSatish Balay 
19437264ac53SSatish Balay   /* if the a->i are the same */
19448108c231SLois Curfman McInnes   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
194577c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
19467264ac53SSatish Balay   }
19477264ac53SSatish Balay 
19487264ac53SSatish Balay   /* if a->j are the same */
1949bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
195077c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1951bcd2baecSBarry Smith   }
1952bcd2baecSBarry Smith 
1953bcd2baecSBarry Smith   /* if a->a are the same */
195419bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
195577c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
19567264ac53SSatish Balay   }
195777c4ece6SBarry Smith   *flg = PETSC_TRUE;
19587264ac53SSatish Balay   return 0;
19597264ac53SSatish Balay 
19607264ac53SSatish Balay }
1961