xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 02ab625ae79bc45e9e37e89987ac6457cedaa9da)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*02ab625aSSatish Balay static char vcid[] = "$Id: aij.c,v 1.235 1997/09/17 23:14:41 bsmith Exp balay $";
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"
158d195f9aSBarry Smith #include "src/inline/dot.h"
16f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1717ab2063SBarry Smith 
18a2ce50c7SBarry Smith /*
19a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
20a2ce50c7SBarry Smith */
215615d1e5SSatish Balay #undef __FUNC__
225615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ"
23a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
24a2ce50c7SBarry Smith {
25a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
26a2ce50c7SBarry Smith   int        ierr = 1;
27a2ce50c7SBarry Smith 
28e3372554SBarry Smith   SETERRQ(ierr,0,"Not implemented");
29a2ce50c7SBarry Smith }
30a2ce50c7SBarry Smith 
31bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3217ab2063SBarry Smith 
335615d1e5SSatish Balay #undef __FUNC__
345615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ"
358f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
366945ee14SBarry Smith                            PetscTruth *done)
3717ab2063SBarry Smith {
38416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
396945ee14SBarry Smith   int        ierr,i,ishift;
4017ab2063SBarry Smith 
4131625ec5SSatish Balay   *m     = A->m;
426945ee14SBarry Smith   if (!ia) return 0;
436945ee14SBarry Smith   ishift = a->indexshift;
446945ee14SBarry Smith   if (symmetric) {
4531625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
466945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
4731625ec5SSatish Balay     int nz = a->i[a->m];
483b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
4931625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
503b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
513b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
5231625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
536945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
5431625ec5SSatish Balay     int nz = a->i[a->m] + 1;
553b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
5631625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
573b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
583b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
5931625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
606945ee14SBarry Smith   } else {
616945ee14SBarry Smith     *ia = a->i; *ja = a->j;
62a2ce50c7SBarry Smith   }
63a2ce50c7SBarry Smith 
64a2744918SBarry Smith   return 0;
65a2744918SBarry Smith }
66a2744918SBarry Smith 
675615d1e5SSatish Balay #undef __FUNC__
685615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ"
698f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
706945ee14SBarry Smith                                PetscTruth *done)
716945ee14SBarry Smith {
726945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
733b2fbd54SBarry Smith   int        ishift = a->indexshift;
746945ee14SBarry Smith 
756945ee14SBarry Smith   if (!ia) return 0;
763b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
776945ee14SBarry Smith     PetscFree(*ia);
786945ee14SBarry Smith     PetscFree(*ja);
79bcd2baecSBarry Smith   }
8017ab2063SBarry Smith   return 0;
8117ab2063SBarry Smith }
8217ab2063SBarry Smith 
835615d1e5SSatish Balay #undef __FUNC__
845615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ"
8543a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
863b2fbd54SBarry Smith                            PetscTruth *done)
873b2fbd54SBarry Smith {
883b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
89a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
90a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
913b2fbd54SBarry Smith 
923b2fbd54SBarry Smith   *nn     = A->n;
933b2fbd54SBarry Smith   if (!ia) return 0;
943b2fbd54SBarry Smith   if (symmetric) {
95179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
963b2fbd54SBarry Smith   } else {
9761d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
983b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
993b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
100a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
1013b2fbd54SBarry Smith     jj = a->j;
1023b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
1033b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
1043b2fbd54SBarry Smith     }
1053b2fbd54SBarry Smith     cia[0] = oshift;
1063b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
1073b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
1083b2fbd54SBarry Smith     }
1093b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1103b2fbd54SBarry Smith     jj = a->j;
111a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
112a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
113a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1143b2fbd54SBarry Smith         col = *jj++ + ishift;
1153b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1163b2fbd54SBarry Smith       }
1173b2fbd54SBarry Smith     }
1183b2fbd54SBarry Smith     PetscFree(collengths);
1193b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1203b2fbd54SBarry Smith   }
1213b2fbd54SBarry Smith 
1223b2fbd54SBarry Smith   return 0;
1233b2fbd54SBarry Smith }
1243b2fbd54SBarry Smith 
1255615d1e5SSatish Balay #undef __FUNC__
1265615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ"
12743a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1283b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1293b2fbd54SBarry Smith {
1303b2fbd54SBarry Smith   if (!ia) return 0;
1313b2fbd54SBarry Smith 
1323b2fbd54SBarry Smith   PetscFree(*ia);
1333b2fbd54SBarry Smith   PetscFree(*ja);
1343b2fbd54SBarry Smith 
1353b2fbd54SBarry Smith   return 0;
1363b2fbd54SBarry Smith }
1373b2fbd54SBarry Smith 
138227d817aSBarry Smith #define CHUNKSIZE   15
13917ab2063SBarry Smith 
1405615d1e5SSatish Balay #undef __FUNC__
1415615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ"
14261d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
14317ab2063SBarry Smith {
144416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
145416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1464b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
147d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
148416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
14917ab2063SBarry Smith 
15017ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
151416022c9SBarry Smith     row  = im[k];
1523b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
153e3372554SBarry Smith     if (row < 0) SETERRQ(1,0,"Negative row");
154e3372554SBarry Smith     if (row >= a->m) SETERRQ(1,0,"Row too large");
1553b2fbd54SBarry Smith #endif
15617ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
15717ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
158416022c9SBarry Smith     low = 0;
15917ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1603b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
161e3372554SBarry Smith       if (in[l] < 0) SETERRQ(1,0,"Negative column");
162e3372554SBarry Smith       if (in[l] >= a->n) SETERRQ(1,0,"Column too large");
1633b2fbd54SBarry Smith #endif
1644b0e389bSBarry Smith       col = in[l] - shift;
1654b0e389bSBarry Smith       if (roworiented) {
1664b0e389bSBarry Smith         value = *v++;
1674b0e389bSBarry Smith       }
1684b0e389bSBarry Smith       else {
1694b0e389bSBarry Smith         value = v[k + l*m];
1704b0e389bSBarry Smith       }
171416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
172416022c9SBarry Smith       while (high-low > 5) {
173416022c9SBarry Smith         t = (low+high)/2;
174416022c9SBarry Smith         if (rp[t] > col) high = t;
175416022c9SBarry Smith         else             low  = t;
17617ab2063SBarry Smith       }
177416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
17817ab2063SBarry Smith         if (rp[i] > col) break;
17917ab2063SBarry Smith         if (rp[i] == col) {
180416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
18117ab2063SBarry Smith           else                  ap[i] = value;
18217ab2063SBarry Smith           goto noinsert;
18317ab2063SBarry Smith         }
18417ab2063SBarry Smith       }
185c2653b3dSLois Curfman McInnes       if (nonew == 1) goto noinsert;
18611ebbc71SLois Curfman McInnes       else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix");
18717ab2063SBarry Smith       if (nrow >= rmax) {
18817ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
189416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
19017ab2063SBarry Smith         Scalar *new_a;
19117ab2063SBarry Smith 
19211ebbc71SLois Curfman McInnes         if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix");
19396854ed6SLois Curfman McInnes 
19417ab2063SBarry Smith         /* malloc new storage space */
195416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1960452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
19717ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
19817ab2063SBarry Smith         new_i   = new_j + new_nz;
19917ab2063SBarry Smith 
20017ab2063SBarry Smith         /* copy over old data into new slots */
20117ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
202416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
203416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
204416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
205416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
20617ab2063SBarry Smith                                                            len*sizeof(int));
207416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
208416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
20917ab2063SBarry Smith                                                            len*sizeof(Scalar));
21017ab2063SBarry Smith         /* free up old matrix storage */
2110452661fSBarry Smith         PetscFree(a->a);
2120452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
213416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
214416022c9SBarry Smith         a->singlemalloc = 1;
21517ab2063SBarry Smith 
21617ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
217416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
218416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
219416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
220b810aeb4SBarry Smith         a->reallocs++;
22117ab2063SBarry Smith       }
222416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
223416022c9SBarry Smith       /* shift up all the later entries in this row */
224416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
22517ab2063SBarry Smith         rp[ii+1] = rp[ii];
22617ab2063SBarry Smith         ap[ii+1] = ap[ii];
22717ab2063SBarry Smith       }
22817ab2063SBarry Smith       rp[i] = col;
22917ab2063SBarry Smith       ap[i] = value;
23017ab2063SBarry Smith       noinsert:;
231416022c9SBarry Smith       low = i + 1;
23217ab2063SBarry Smith     }
23317ab2063SBarry Smith     ailen[row] = nrow;
23417ab2063SBarry Smith   }
23517ab2063SBarry Smith   return 0;
23617ab2063SBarry Smith }
23717ab2063SBarry Smith 
2385615d1e5SSatish Balay #undef __FUNC__
2395615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ"
2408f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2417eb43aa7SLois Curfman McInnes {
2427eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
243b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2447eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2457eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2467eb43aa7SLois Curfman McInnes 
2477eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2487eb43aa7SLois Curfman McInnes     row  = im[k];
249e3372554SBarry Smith     if (row < 0) SETERRQ(1,0,"Negative row");
250e3372554SBarry Smith     if (row >= a->m) SETERRQ(1,0,"Row too large");
2517eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2527eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2537eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
254e3372554SBarry Smith       if (in[l] < 0) SETERRQ(1,0,"Negative column");
255e3372554SBarry Smith       if (in[l] >= a->n) SETERRQ(1,0,"Column too large");
2567eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2577eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2587eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2597eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2607eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2617eb43aa7SLois Curfman McInnes         else             low  = t;
2627eb43aa7SLois Curfman McInnes       }
2637eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2647eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2657eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
266b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2677eb43aa7SLois Curfman McInnes           goto finished;
2687eb43aa7SLois Curfman McInnes         }
2697eb43aa7SLois Curfman McInnes       }
270b49de8d1SLois Curfman McInnes       *v++ = zero;
2717eb43aa7SLois Curfman McInnes       finished:;
2727eb43aa7SLois Curfman McInnes     }
2737eb43aa7SLois Curfman McInnes   }
2747eb43aa7SLois Curfman McInnes   return 0;
2757eb43aa7SLois Curfman McInnes }
2767eb43aa7SLois Curfman McInnes 
27717ab2063SBarry Smith 
2785615d1e5SSatish Balay #undef __FUNC__
2795615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary"
2808f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
28117ab2063SBarry Smith {
282416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
283416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
28417ab2063SBarry Smith 
28590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2860452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
287416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
288416022c9SBarry Smith   col_lens[1] = a->m;
289416022c9SBarry Smith   col_lens[2] = a->n;
290416022c9SBarry Smith   col_lens[3] = a->nz;
291416022c9SBarry Smith 
292416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
293416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
294416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
29517ab2063SBarry Smith   }
29677c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
2970452661fSBarry Smith   PetscFree(col_lens);
298416022c9SBarry Smith 
299416022c9SBarry Smith   /* store column indices (zero start index) */
300416022c9SBarry Smith   if (a->indexshift) {
301416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
30217ab2063SBarry Smith   }
30377c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr);
304416022c9SBarry Smith   if (a->indexshift) {
305416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
30617ab2063SBarry Smith   }
307416022c9SBarry Smith 
308416022c9SBarry Smith   /* store nonzero values */
30977c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
31017ab2063SBarry Smith   return 0;
31117ab2063SBarry Smith }
312416022c9SBarry Smith 
3135615d1e5SSatish Balay #undef __FUNC__
3145615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII"
3158f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
316416022c9SBarry Smith {
317416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
318496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
31917ab2063SBarry Smith   FILE        *fd;
32017ab2063SBarry Smith   char        *outputname;
32117ab2063SBarry Smith 
32290ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
323416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
32490ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
325a93ec695SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
32695e01e2fSLois Curfman McInnes     return 0;
32795e01e2fSLois Curfman McInnes   }
328a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
329496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
330496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
331496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
33295e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
33395e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
33417ab2063SBarry Smith   }
335a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
336d00d2cf4SBarry Smith     int nofinalvalue = 0;
337d00d2cf4SBarry Smith     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) {
338d00d2cf4SBarry Smith       nofinalvalue = 1;
339d00d2cf4SBarry Smith     }
340416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3414e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
342d00d2cf4SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue);
34317ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
34417ab2063SBarry Smith 
34517ab2063SBarry Smith     for (i=0; i<m; i++) {
346416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
34717ab2063SBarry Smith #if defined(PETSC_COMPLEX)
3486945ee14SBarry Smith         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]),
349416022c9SBarry Smith                    imag(a->a[j]));
35017ab2063SBarry Smith #else
3517a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
35217ab2063SBarry Smith #endif
35317ab2063SBarry Smith       }
35417ab2063SBarry Smith     }
355d00d2cf4SBarry Smith     if (nofinalvalue) {
356d00d2cf4SBarry Smith       fprintf(fd,"%d %d  %18.16e\n", m, a->n, 0.0);
357d00d2cf4SBarry Smith     }
35817ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
35917ab2063SBarry Smith   }
360a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
36144cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
36244cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
36344cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
36444cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
365766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0 && real(a->a[j]) != 0.0)
36644cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
367766eeae4SLois Curfman McInnes         else if (imag(a->a[j]) < 0.0 && real(a->a[j]) != 0.0)
368766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
36944cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
37044cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
37144cd7ae7SLois Curfman McInnes #else
37244cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
37344cd7ae7SLois Curfman McInnes #endif
37444cd7ae7SLois Curfman McInnes       }
37544cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
37644cd7ae7SLois Curfman McInnes     }
37744cd7ae7SLois Curfman McInnes   }
378496be53dSLois Curfman McInnes   else if (format == VIEWER_FORMAT_ASCII_SYMMODU) {
379496be53dSLois Curfman McInnes     int nzd=0, fshift=1, *sptr;
3802e44a96cSLois Curfman McInnes     sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr);
381496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
382496be53dSLois Curfman McInnes       sptr[i] = nzd+1;
383496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
384496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
385496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX)
386496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) nzd++;
387496be53dSLois Curfman McInnes #else
388496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) nzd++;
389496be53dSLois Curfman McInnes #endif
390496be53dSLois Curfman McInnes         }
391496be53dSLois Curfman McInnes       }
392496be53dSLois Curfman McInnes     }
3932e44a96cSLois Curfman McInnes     sptr[m] = nzd+1;
394496be53dSLois Curfman McInnes     fprintf(fd," %d %d\n\n",m,nzd);
3952e44a96cSLois Curfman McInnes     for ( i=0; i<m+1; i+=6 ) {
3962e44a96cSLois 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]);
3972e44a96cSLois 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]);
3982e44a96cSLois 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]);
3992e44a96cSLois Curfman McInnes       else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]);
4002e44a96cSLois Curfman McInnes       else if (i<m)   fprintf(fd," %d %d\n",sptr[i],sptr[i+1]);
4017272d637SLois Curfman McInnes       else            fprintf(fd," %d\n",sptr[i]);
402496be53dSLois Curfman McInnes     }
403496be53dSLois Curfman McInnes     fprintf(fd,"\n");
404496be53dSLois Curfman McInnes     PetscFree(sptr);
405496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
406496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
407496be53dSLois Curfman McInnes         if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift);
408496be53dSLois Curfman McInnes       }
409496be53dSLois Curfman McInnes       fprintf(fd,"\n");
410496be53dSLois Curfman McInnes     }
411496be53dSLois Curfman McInnes     fprintf(fd,"\n");
412496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
413496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
414496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
415496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX)
416496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0)
417496be53dSLois Curfman McInnes             fprintf(fd," %18.16e %18.16e ",real(a->a[j]),imag(a->a[j]));
418496be53dSLois Curfman McInnes #else
419496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]);
420496be53dSLois Curfman McInnes #endif
421496be53dSLois Curfman McInnes         }
422496be53dSLois Curfman McInnes       }
423496be53dSLois Curfman McInnes       fprintf(fd,"\n");
424496be53dSLois Curfman McInnes     }
425496be53dSLois Curfman McInnes   }
42617ab2063SBarry Smith   else {
42717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
42817ab2063SBarry Smith       fprintf(fd,"row %d:",i);
429416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
43017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
431766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0) {
432416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
433766eeae4SLois Curfman McInnes         } else if (imag(a->a[j]) < 0.0) {
434766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
43517ab2063SBarry Smith         }
43617ab2063SBarry Smith         else {
437416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
43817ab2063SBarry Smith         }
43917ab2063SBarry Smith #else
440416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
44117ab2063SBarry Smith #endif
44217ab2063SBarry Smith       }
44317ab2063SBarry Smith       fprintf(fd,"\n");
44417ab2063SBarry Smith     }
44517ab2063SBarry Smith   }
44617ab2063SBarry Smith   fflush(fd);
447416022c9SBarry Smith   return 0;
448416022c9SBarry Smith }
449416022c9SBarry Smith 
4505615d1e5SSatish Balay #undef __FUNC__
4515615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Draw"
4528f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
453416022c9SBarry Smith {
454416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
455cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
4560513a670SBarry Smith   int         format;
45794a9d846SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r,maxv = 0.0;
458bcd2baecSBarry Smith   Draw        draw;
459cddf8d76SBarry Smith   DrawButton  button;
46019bcc07fSBarry Smith   PetscTruth  isnull;
461cddf8d76SBarry Smith 
4620513a670SBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
4630513a670SBarry Smith   ierr = DrawClear(draw); CHKERRQ(ierr);
4640513a670SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
46519bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
46619bcc07fSBarry Smith 
467416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
468416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
469416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
470416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
4710513a670SBarry Smith 
4720513a670SBarry Smith   if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4730513a670SBarry Smith     /* Blue for negative, Cyan for zero and  Red for positive */
474cddf8d76SBarry Smith     color = DRAW_BLUE;
475416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
476cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
477416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
478cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
479cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
480cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
481cddf8d76SBarry Smith #else
482cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
483cddf8d76SBarry Smith #endif
484cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
485cddf8d76SBarry Smith       }
486cddf8d76SBarry Smith     }
487cddf8d76SBarry Smith     color = DRAW_CYAN;
488cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
489cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
490cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
491cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
492cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
493cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
494cddf8d76SBarry Smith       }
495cddf8d76SBarry Smith     }
496cddf8d76SBarry Smith     color = DRAW_RED;
497cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
498cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
499cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
500cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
501cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
502cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
503cddf8d76SBarry Smith #else
504cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
505cddf8d76SBarry Smith #endif
506cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
507416022c9SBarry Smith       }
508416022c9SBarry Smith     }
5090513a670SBarry Smith   } else {
5100513a670SBarry Smith     /* use contour shading to indicate magnitude of values */
5110513a670SBarry Smith     /* first determine max of all nonzero values */
5120513a670SBarry Smith     int    nz = a->nz,count;
5130513a670SBarry Smith     Draw   popup;
5140513a670SBarry Smith 
5150513a670SBarry Smith     for ( i=0; i<nz; i++ ) {
5160513a670SBarry Smith       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
5170513a670SBarry Smith     }
5180513a670SBarry Smith     ierr = DrawCreatePopUp(draw,&popup); CHKERRQ(ierr);
5190513a670SBarry Smith     ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr);
5200513a670SBarry Smith     count = 0;
5210513a670SBarry Smith     for ( i=0; i<m; i++ ) {
5220513a670SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
5230513a670SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5240513a670SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5250513a670SBarry Smith         color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5260513a670SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5270513a670SBarry Smith         count++;
5280513a670SBarry Smith       }
5290513a670SBarry Smith     }
5300513a670SBarry Smith   }
531416022c9SBarry Smith   DrawFlush(draw);
532cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
533cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
534cddf8d76SBarry Smith 
535cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
5366945ee14SBarry Smith   ierr = DrawCheckResizedWindow(draw);
537cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
538cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
539cddf8d76SBarry Smith     DrawClear(draw);
540cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
541cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
542cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
543cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
544cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
545cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
546cddf8d76SBarry Smith     w *= scale; h *= scale;
547cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5480513a670SBarry Smith     if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
5490513a670SBarry Smith       /* Blue for negative, Cyan for zero and  Red for positive */
550cddf8d76SBarry Smith       color = DRAW_BLUE;
551cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
552cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
553cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
554cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
555cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
556cddf8d76SBarry Smith           if (real(a->a[j]) >=  0.) continue;
557cddf8d76SBarry Smith #else
558cddf8d76SBarry Smith           if (a->a[j] >=  0.) continue;
559cddf8d76SBarry Smith #endif
560cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
561cddf8d76SBarry Smith         }
562cddf8d76SBarry Smith       }
563cddf8d76SBarry Smith       color = DRAW_CYAN;
564cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
565cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
566cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
567cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
568cddf8d76SBarry Smith           if (a->a[j] !=  0.) continue;
569cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
570cddf8d76SBarry Smith         }
571cddf8d76SBarry Smith       }
572cddf8d76SBarry Smith       color = DRAW_RED;
573cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
574cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
575cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
576cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
577cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
578cddf8d76SBarry Smith           if (real(a->a[j]) <=  0.) continue;
579cddf8d76SBarry Smith #else
580cddf8d76SBarry Smith           if (a->a[j] <=  0.) continue;
581cddf8d76SBarry Smith #endif
582cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
583cddf8d76SBarry Smith         }
584cddf8d76SBarry Smith       }
5850513a670SBarry Smith     } else {
5860513a670SBarry Smith       /* use contour shading to indicate magnitude of values */
5870513a670SBarry Smith       int count = 0;
5880513a670SBarry Smith       for ( i=0; i<m; i++ ) {
5890513a670SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
5900513a670SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5910513a670SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
5920513a670SBarry Smith           color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5930513a670SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5940513a670SBarry Smith           count++;
5950513a670SBarry Smith         }
5960513a670SBarry Smith       }
5970513a670SBarry Smith     }
5980513a670SBarry Smith 
5996945ee14SBarry Smith     ierr = DrawCheckResizedWindow(draw);
600cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
601cddf8d76SBarry Smith   }
602416022c9SBarry Smith   return 0;
603416022c9SBarry Smith }
604416022c9SBarry Smith 
6055615d1e5SSatish Balay #undef __FUNC__
606d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ"
6078f6be9afSLois Curfman McInnes int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
608416022c9SBarry Smith {
609416022c9SBarry Smith   Mat         A = (Mat) obj;
610416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
611bcd2baecSBarry Smith   ViewerType  vtype;
612bcd2baecSBarry Smith   int         ierr;
613416022c9SBarry Smith 
614bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
615bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
616416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
617416022c9SBarry Smith   }
618bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
619416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
620416022c9SBarry Smith   }
621bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
622416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
623416022c9SBarry Smith   }
624bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
625bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
62617ab2063SBarry Smith   }
62717ab2063SBarry Smith   return 0;
62817ab2063SBarry Smith }
62919bcc07fSBarry Smith 
630c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
6315615d1e5SSatish Balay #undef __FUNC__
6325615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ"
6338f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
63417ab2063SBarry Smith {
635416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
63641c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
63743ee02c3SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0;
638416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
63917ab2063SBarry Smith 
6406d4a8577SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
64117ab2063SBarry Smith 
64243ee02c3SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
64317ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
644416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
64517ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
64694a9d846SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
64717ab2063SBarry Smith     if (fshift) {
648416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
64917ab2063SBarry Smith       N = ailen[i];
65017ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
65117ab2063SBarry Smith         ip[j-fshift] = ip[j];
65217ab2063SBarry Smith         ap[j-fshift] = ap[j];
65317ab2063SBarry Smith       }
65417ab2063SBarry Smith     }
65517ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
65617ab2063SBarry Smith   }
65717ab2063SBarry Smith   if (m) {
65817ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
65917ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
66017ab2063SBarry Smith   }
66117ab2063SBarry Smith   /* reset ilen and imax for each row */
66217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
66317ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
66417ab2063SBarry Smith   }
665416022c9SBarry Smith   a->nz = ai[m] + shift;
66617ab2063SBarry Smith 
66717ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
668416022c9SBarry Smith   if (fshift && a->diag) {
6690452661fSBarry Smith     PetscFree(a->diag);
670416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
671416022c9SBarry Smith     a->diag = 0;
67217ab2063SBarry Smith   }
6734e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
6744e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
6754e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
676b810aeb4SBarry Smith            a->reallocs);
67794a9d846SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax);
678dd5f02e7SSatish Balay   a->reallocs          = 0;
6794e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
6804e220ebcSLois Curfman McInnes 
68176dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
68241c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
68317ab2063SBarry Smith   return 0;
68417ab2063SBarry Smith }
68517ab2063SBarry Smith 
6865615d1e5SSatish Balay #undef __FUNC__
6875615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ"
6888f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A)
68917ab2063SBarry Smith {
690416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
691cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
69217ab2063SBarry Smith   return 0;
69317ab2063SBarry Smith }
694416022c9SBarry Smith 
6955615d1e5SSatish Balay #undef __FUNC__
6965615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ"
69717ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
69817ab2063SBarry Smith {
699416022c9SBarry Smith   Mat        A  = (Mat) obj;
700416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
701d5d45c9bSBarry Smith 
70217ab2063SBarry Smith #if defined(PETSC_LOG)
703416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
70417ab2063SBarry Smith #endif
7050452661fSBarry Smith   PetscFree(a->a);
7060452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
7070452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
7080452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
7090452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
7100452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
71176dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
7120452661fSBarry Smith   PetscFree(a);
713eed86810SBarry Smith 
714f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
715f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
71617ab2063SBarry Smith   return 0;
71717ab2063SBarry Smith }
71817ab2063SBarry Smith 
7195615d1e5SSatish Balay #undef __FUNC__
7205615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
7218f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A)
72217ab2063SBarry Smith {
72317ab2063SBarry Smith   return 0;
72417ab2063SBarry Smith }
72517ab2063SBarry Smith 
7265615d1e5SSatish Balay #undef __FUNC__
7275615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
7288f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op)
72917ab2063SBarry Smith {
730416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7316d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
7326d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
7336d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
734219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
7356d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
736c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
73796854ed6SLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
7386d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
7396d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
740219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7416d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7426d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
74390f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
7442b362799SSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES)
74594a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
7466d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
747e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");}
7486d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
7496d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
7506d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
7516d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
7526d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
753e2f28af5SBarry Smith   else
754e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
75517ab2063SBarry Smith   return 0;
75617ab2063SBarry Smith }
75717ab2063SBarry Smith 
7585615d1e5SSatish Balay #undef __FUNC__
7595615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
7608f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
76117ab2063SBarry Smith {
762416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
763416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
76417ab2063SBarry Smith   Scalar     *x, zero = 0.0;
76517ab2063SBarry Smith 
76617ab2063SBarry Smith   VecSet(&zero,v);
76790f02eecSBarry Smith   VecGetArray_Fast(v,x); VecGetLocalSize(v,&n);
768e3372554SBarry Smith   if (n != a->m) SETERRQ(1,0,"Nonconforming matrix and vector");
769416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
770416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
771416022c9SBarry Smith       if (a->j[j]+shift == i) {
772416022c9SBarry Smith         x[i] = a->a[j];
77317ab2063SBarry Smith         break;
77417ab2063SBarry Smith       }
77517ab2063SBarry Smith     }
77617ab2063SBarry Smith   }
77717ab2063SBarry Smith   return 0;
77817ab2063SBarry Smith }
77917ab2063SBarry Smith 
78017ab2063SBarry Smith /* -------------------------------------------------------*/
78117ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
78217ab2063SBarry Smith /* -------------------------------------------------------*/
7835615d1e5SSatish Balay #undef __FUNC__
7845615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
78544cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
78617ab2063SBarry Smith {
787416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
78817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
789416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
79017ab2063SBarry Smith 
79190f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
792cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
79317ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
79417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
795416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
796416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
797416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
79817ab2063SBarry Smith     alpha = x[i];
79917ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
80017ab2063SBarry Smith   }
801416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
80217ab2063SBarry Smith   return 0;
80317ab2063SBarry Smith }
804d5d45c9bSBarry Smith 
8055615d1e5SSatish Balay #undef __FUNC__
8065615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
80744cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
80817ab2063SBarry Smith {
809416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
81017ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
811416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
81217ab2063SBarry Smith 
81390f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
81417ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
81517ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
81617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
817416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
818416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
819416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
82017ab2063SBarry Smith     alpha = x[i];
82117ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
82217ab2063SBarry Smith   }
82390f02eecSBarry Smith   PLogFlops(2*a->nz);
82417ab2063SBarry Smith   return 0;
82517ab2063SBarry Smith }
82617ab2063SBarry Smith 
8275615d1e5SSatish Balay #undef __FUNC__
8285615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
82944cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
83017ab2063SBarry Smith {
831416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
83217ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
8339ea0dfa2SSatish Balay   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j;
83417ab2063SBarry Smith 
83590f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
83617ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
837416022c9SBarry Smith   idx  = a->j;
8389b88f4a7SBarry Smith   v    = a->a + shift; /* shift for Fortran start by 1 indexing */
839416022c9SBarry Smith   ii   = a->i;
8408d195f9aSBarry Smith #if defined(USE_FORTRAN_KERNELS)
8418d195f9aSBarry Smith   fortranmultaij_(&m,x,ii,idx,v,y);
8428d195f9aSBarry Smith #else
84317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8449ea0dfa2SSatish Balay     jrow = ii[i];
8459ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
84617ab2063SBarry Smith     sum  = 0.0;
8479ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8489ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8499ea0dfa2SSatish Balay      }
85017ab2063SBarry Smith     y[i] = sum;
85117ab2063SBarry Smith   }
8528d195f9aSBarry Smith #endif
853416022c9SBarry Smith   PLogFlops(2*a->nz - m);
85417ab2063SBarry Smith   return 0;
85517ab2063SBarry Smith }
85617ab2063SBarry Smith 
8575615d1e5SSatish Balay #undef __FUNC__
8585615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
85944cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
86017ab2063SBarry Smith {
861416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
86217ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
8639ea0dfa2SSatish Balay   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j;
8649ea0dfa2SSatish Balay 
86517ab2063SBarry Smith 
86690f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); VecGetArray_Fast(zz,z);
86717ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
868cddf8d76SBarry Smith   idx  = a->j;
8699b88f4a7SBarry Smith   v    = a->a + shift; /* shift for Fortran start by 1 indexing */
870cddf8d76SBarry Smith   ii   = a->i;
871*02ab625aSSatish Balay #if defined(USE_FORTRAN_KERNELS)
872*02ab625aSSatish Balay   fortranmultaddaij_(&m,x,ii,idx,v,y,z);
873*02ab625aSSatish Balay #else
87417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8759ea0dfa2SSatish Balay     jrow = ii[i];
8769ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
87717ab2063SBarry Smith     sum  = y[i];
8789ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8799ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8809ea0dfa2SSatish Balay      }
88117ab2063SBarry Smith     z[i] = sum;
88217ab2063SBarry Smith   }
883*02ab625aSSatish Balay #endif
884416022c9SBarry Smith   PLogFlops(2*a->nz);
88517ab2063SBarry Smith   return 0;
88617ab2063SBarry Smith }
88717ab2063SBarry Smith 
88817ab2063SBarry Smith /*
88917ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
89017ab2063SBarry Smith */
89117ab2063SBarry Smith 
8925615d1e5SSatish Balay #undef __FUNC__
8935615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
894416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
89517ab2063SBarry Smith {
896416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
897416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
89817ab2063SBarry Smith 
8990452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
900416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
901416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
902416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
903416022c9SBarry Smith       if (a->j[j]+shift == i) {
90417ab2063SBarry Smith         diag[i] = j - shift;
90517ab2063SBarry Smith         break;
90617ab2063SBarry Smith       }
90717ab2063SBarry Smith     }
90817ab2063SBarry Smith   }
909416022c9SBarry Smith   a->diag = diag;
91017ab2063SBarry Smith   return 0;
91117ab2063SBarry Smith }
91217ab2063SBarry Smith 
9135615d1e5SSatish Balay #undef __FUNC__
9145615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
91544cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
91617ab2063SBarry Smith                            double fshift,int its,Vec xx)
91717ab2063SBarry Smith {
918416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
919416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
920d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
92117ab2063SBarry Smith 
92290f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(bb,b);
923d00d2cf4SBarry Smith   if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);}
924416022c9SBarry Smith   diag = a->diag;
925416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
92617ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
92717ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
92817ab2063SBarry Smith     bs = b + shift;
92917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
930416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
931416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
932416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
933416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
93417ab2063SBarry Smith         sum  = b[i]*d/omega;
93517ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
93617ab2063SBarry Smith         x[i] = sum;
93717ab2063SBarry Smith     }
93817ab2063SBarry Smith     return 0;
93917ab2063SBarry Smith   }
94017ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
941e3372554SBarry Smith     SETERRQ(1,0,"SOR_APPLY_LOWER is not done");
94217ab2063SBarry Smith   }
943416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
94417ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
94517ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
94617ab2063SBarry Smith 
94717ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
94817ab2063SBarry Smith 
94917ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
95017ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
95117ab2063SBarry Smith     is the relaxation factor.
95217ab2063SBarry Smith     */
9530452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
95417ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
95517ab2063SBarry Smith 
95617ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
95717ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
958416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
959416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
960416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
961416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
96217ab2063SBarry Smith       sum  = b[i];
96317ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
96417ab2063SBarry Smith       x[i] = omega*(sum/d);
96517ab2063SBarry Smith     }
96617ab2063SBarry Smith 
96717ab2063SBarry Smith     /*  t = b - (2*E - D)x */
968416022c9SBarry Smith     v = a->a;
96917ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
97017ab2063SBarry Smith 
97117ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
972416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
973416022c9SBarry Smith     diag = a->diag;
97417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
975416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
976416022c9SBarry Smith       n    = diag[i] - a->i[i];
977416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
978416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
97917ab2063SBarry Smith       sum  = t[i];
98017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
98117ab2063SBarry Smith       t[i] = omega*(sum/d);
98217ab2063SBarry Smith     }
98317ab2063SBarry Smith 
98417ab2063SBarry Smith     /*  x = x + t */
98517ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
9860452661fSBarry Smith     PetscFree(t);
98717ab2063SBarry Smith     return 0;
98817ab2063SBarry Smith   }
98917ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
99017ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
99117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
992416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
993416022c9SBarry Smith         n    = diag[i] - a->i[i];
994416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
995416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
99617ab2063SBarry Smith         sum  = b[i];
99717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
99817ab2063SBarry Smith         x[i] = omega*(sum/d);
99917ab2063SBarry Smith       }
100017ab2063SBarry Smith       xb = x;
100117ab2063SBarry Smith     }
100217ab2063SBarry Smith     else xb = b;
100317ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
100417ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
100517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1006416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
100717ab2063SBarry Smith       }
100817ab2063SBarry Smith     }
100917ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
101017ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1011416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1012416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
1013416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1014416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
101517ab2063SBarry Smith         sum  = xb[i];
101617ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
101717ab2063SBarry Smith         x[i] = omega*(sum/d);
101817ab2063SBarry Smith       }
101917ab2063SBarry Smith     }
102017ab2063SBarry Smith     its--;
102117ab2063SBarry Smith   }
102217ab2063SBarry Smith   while (its--) {
102317ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
102417ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1025416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1026416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1027416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1028416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
102917ab2063SBarry Smith         sum  = b[i];
103017ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10317e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
103217ab2063SBarry Smith       }
103317ab2063SBarry Smith     }
103417ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
103517ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1036416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1037416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1038416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1039416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
104017ab2063SBarry Smith         sum  = b[i];
104117ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10427e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
104317ab2063SBarry Smith       }
104417ab2063SBarry Smith     }
104517ab2063SBarry Smith   }
104617ab2063SBarry Smith   return 0;
104717ab2063SBarry Smith }
104817ab2063SBarry Smith 
10495615d1e5SSatish Balay #undef __FUNC__
10505615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10518f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
105217ab2063SBarry Smith {
1053416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
10544e220ebcSLois Curfman McInnes 
10554e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
10564e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
10574e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
10584e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
10594e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10604e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
10614e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
10624e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
10634e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
10644e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
10654e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
10664e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
10674e220ebcSLois Curfman McInnes   info->memory         = A->mem;
10684e220ebcSLois Curfman McInnes   if (A->factor) {
10694e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10704e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
10714e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
10724e220ebcSLois Curfman McInnes   } else {
10734e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
10744e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
10754e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
10764e220ebcSLois Curfman McInnes   }
107717ab2063SBarry Smith   return 0;
107817ab2063SBarry Smith }
107917ab2063SBarry Smith 
108017ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
108117ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
108217ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
108317ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
108417ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
108517ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
108617ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
108717ab2063SBarry Smith 
10885615d1e5SSatish Balay #undef __FUNC__
10895615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
10908f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
109117ab2063SBarry Smith {
1092416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1093416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
109417ab2063SBarry Smith 
109577c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
109617ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
109717ab2063SBarry Smith   if (diag) {
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       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1101416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1102416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1103416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
110417ab2063SBarry Smith       }
110517ab2063SBarry Smith       else {
110617ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
110717ab2063SBarry Smith         CHKERRQ(ierr);
110817ab2063SBarry Smith       }
110917ab2063SBarry Smith     }
111017ab2063SBarry Smith   }
111117ab2063SBarry Smith   else {
111217ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1113e3372554SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range");
1114416022c9SBarry Smith       a->ilen[rows[i]] = 0;
111517ab2063SBarry Smith     }
111617ab2063SBarry Smith   }
111717ab2063SBarry Smith   ISRestoreIndices(is,&rows);
111843a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
111917ab2063SBarry Smith   return 0;
112017ab2063SBarry Smith }
112117ab2063SBarry Smith 
11225615d1e5SSatish Balay #undef __FUNC__
11235615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11248f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
112517ab2063SBarry Smith {
1126416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1127416022c9SBarry Smith   *m = a->m; *n = a->n;
112817ab2063SBarry Smith   return 0;
112917ab2063SBarry Smith }
113017ab2063SBarry Smith 
11315615d1e5SSatish Balay #undef __FUNC__
11325615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11338f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
113417ab2063SBarry Smith {
1135416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1136416022c9SBarry Smith   *m = 0; *n = a->m;
113717ab2063SBarry Smith   return 0;
113817ab2063SBarry Smith }
1139026e39d0SSatish Balay 
11405615d1e5SSatish Balay #undef __FUNC__
11415615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11424e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
114317ab2063SBarry Smith {
1144416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1145c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
114617ab2063SBarry Smith 
1147e3372554SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,0,"Row out of range");
114817ab2063SBarry Smith 
1149416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1150416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
115117ab2063SBarry Smith   if (idx) {
1152416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
11534e093b46SBarry Smith     if (*nz && shift) {
11540452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
115517ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
11564e093b46SBarry Smith     } else if (*nz) {
11574e093b46SBarry Smith       *idx = itmp;
115817ab2063SBarry Smith     }
115917ab2063SBarry Smith     else *idx = 0;
116017ab2063SBarry Smith   }
116117ab2063SBarry Smith   return 0;
116217ab2063SBarry Smith }
116317ab2063SBarry Smith 
11645615d1e5SSatish Balay #undef __FUNC__
11655615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
11664e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
116717ab2063SBarry Smith {
11684e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11694e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
117017ab2063SBarry Smith   return 0;
117117ab2063SBarry Smith }
117217ab2063SBarry Smith 
11735615d1e5SSatish Balay #undef __FUNC__
11745615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
11758f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
117617ab2063SBarry Smith {
1177416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1178416022c9SBarry Smith   Scalar     *v = a->a;
117917ab2063SBarry Smith   double     sum = 0.0;
1180416022c9SBarry Smith   int        i, j,shift = a->indexshift;
118117ab2063SBarry Smith 
118217ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1183416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
118417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
118517ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
118617ab2063SBarry Smith #else
118717ab2063SBarry Smith       sum += (*v)*(*v); v++;
118817ab2063SBarry Smith #endif
118917ab2063SBarry Smith     }
119017ab2063SBarry Smith     *norm = sqrt(sum);
119117ab2063SBarry Smith   }
119217ab2063SBarry Smith   else if (type == NORM_1) {
119317ab2063SBarry Smith     double *tmp;
1194416022c9SBarry Smith     int    *jj = a->j;
119566963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1196cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
119717ab2063SBarry Smith     *norm = 0.0;
1198416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1199a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
120017ab2063SBarry Smith     }
1201416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
120217ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
120317ab2063SBarry Smith     }
12040452661fSBarry Smith     PetscFree(tmp);
120517ab2063SBarry Smith   }
120617ab2063SBarry Smith   else if (type == NORM_INFINITY) {
120717ab2063SBarry Smith     *norm = 0.0;
1208416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1209416022c9SBarry Smith       v = a->a + a->i[j] + shift;
121017ab2063SBarry Smith       sum = 0.0;
1211416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1212cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
121317ab2063SBarry Smith       }
121417ab2063SBarry Smith       if (sum > *norm) *norm = sum;
121517ab2063SBarry Smith     }
121617ab2063SBarry Smith   }
121717ab2063SBarry Smith   else {
1218e3372554SBarry Smith     SETERRQ(1,0,"No support for two norm yet");
121917ab2063SBarry Smith   }
122017ab2063SBarry Smith   return 0;
122117ab2063SBarry Smith }
122217ab2063SBarry Smith 
12235615d1e5SSatish Balay #undef __FUNC__
12245615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12258f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
122617ab2063SBarry Smith {
1227416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1228416022c9SBarry Smith   Mat        C;
1229416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1230416022c9SBarry Smith   int        shift = a->indexshift;
1231d5d45c9bSBarry Smith   Scalar     *array = a->a;
123217ab2063SBarry Smith 
12333638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
1234e3372554SBarry Smith     SETERRQ(1,0,"Square matrix only for in-place");
12350452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1236cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
123717ab2063SBarry Smith   if (shift) {
123817ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
123917ab2063SBarry Smith   }
124017ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1241416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12420452661fSBarry Smith   PetscFree(col);
124317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
124417ab2063SBarry Smith     len = ai[i+1]-ai[i];
1245416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
124617ab2063SBarry Smith     array += len; aj += len;
124717ab2063SBarry Smith   }
124817ab2063SBarry Smith   if (shift) {
124917ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
125017ab2063SBarry Smith   }
125117ab2063SBarry Smith 
12526d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12536d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
125417ab2063SBarry Smith 
12553638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1256416022c9SBarry Smith     *B = C;
125717ab2063SBarry Smith   } else {
1258416022c9SBarry Smith     /* This isn't really an in-place transpose */
12590452661fSBarry Smith     PetscFree(a->a);
12600452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
12610452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
12620452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
12630452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
12640452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
12651073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
12660452661fSBarry Smith     PetscFree(a);
1267f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
12680452661fSBarry Smith     PetscHeaderDestroy(C);
126917ab2063SBarry Smith   }
127017ab2063SBarry Smith   return 0;
127117ab2063SBarry Smith }
127217ab2063SBarry Smith 
12735615d1e5SSatish Balay #undef __FUNC__
12745615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
12758f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
127617ab2063SBarry Smith {
1277416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
127817ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1279d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
128017ab2063SBarry Smith 
128117ab2063SBarry Smith   if (ll) {
12823ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
12833ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
12849b1297e1SSatish Balay     VecGetLocalSize_Fast(ll,m);
1285e3372554SBarry Smith     if (m != a->m) SETERRQ(1,0,"Left scaling vector wrong length");
128690f02eecSBarry Smith     VecGetArray_Fast(ll,l);
1287416022c9SBarry Smith     v = a->a;
128817ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
128917ab2063SBarry Smith       x = l[i];
1290416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
129117ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
129217ab2063SBarry Smith     }
129344cd7ae7SLois Curfman McInnes     PLogFlops(nz);
129417ab2063SBarry Smith   }
129517ab2063SBarry Smith   if (rr) {
12969b1297e1SSatish Balay     VecGetLocalSize_Fast(rr,n);
1297e3372554SBarry Smith     if (n != a->n) SETERRQ(1,0,"Right scaling vector wrong length");
129890f02eecSBarry Smith     VecGetArray_Fast(rr,r);
1299416022c9SBarry Smith     v = a->a; jj = a->j;
130017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
130117ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
130217ab2063SBarry Smith     }
130344cd7ae7SLois Curfman McInnes     PLogFlops(nz);
130417ab2063SBarry Smith   }
130517ab2063SBarry Smith   return 0;
130617ab2063SBarry Smith }
130717ab2063SBarry Smith 
13085615d1e5SSatish Balay #undef __FUNC__
13095615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
13108f6be9afSLois Curfman McInnes int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
131117ab2063SBarry Smith {
1312db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
131302834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
131499141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1315a2744918SBarry Smith   register int sum,lensi;
131699141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
131799141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
131899141d43SSatish Balay   Scalar       *a_new,*mat_a;
1319416022c9SBarry Smith   Mat          C;
132017ab2063SBarry Smith 
1321b48a1e75SSatish Balay   ierr = ISSorted(isrow,(PetscTruth*)&i);
1322e3372554SBarry Smith   if (!i) SETERRQ(1,0,"ISrow is not sorted");
132399141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
1324e3372554SBarry Smith   if (!i) SETERRQ(1,0,"IScol is not sorted");
132599141d43SSatish Balay 
132617ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
132717ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
132817ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
132917ab2063SBarry Smith 
13307264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
133102834360SBarry Smith     /* special case of contiguous rows */
133257faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
133302834360SBarry Smith     starts = lens + ncols;
133402834360SBarry Smith     /* loop over new rows determining lens and starting points */
133502834360SBarry Smith     for (i=0; i<nrows; i++) {
1336a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1337a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
133802834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1339d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
134002834360SBarry Smith           starts[i] = k;
134102834360SBarry Smith           break;
134202834360SBarry Smith 	}
134302834360SBarry Smith       }
1344a2744918SBarry Smith       sum = 0;
134502834360SBarry Smith       while (k < kend) {
1346d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1347a2744918SBarry Smith         sum++;
134802834360SBarry Smith       }
1349a2744918SBarry Smith       lens[i] = sum;
135002834360SBarry Smith     }
135102834360SBarry Smith     /* create submatrix */
1352cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
135308480c60SBarry Smith       int n_cols,n_rows;
135408480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1355e3372554SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,0,"");
1356d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
135708480c60SBarry Smith       C = *B;
135808480c60SBarry Smith     }
135908480c60SBarry Smith     else {
136002834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
136108480c60SBarry Smith     }
1362db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1363db02288aSLois Curfman McInnes 
136402834360SBarry Smith     /* loop over rows inserting into submatrix */
1365db02288aSLois Curfman McInnes     a_new    = c->a;
1366db02288aSLois Curfman McInnes     j_new    = c->j;
1367db02288aSLois Curfman McInnes     i_new    = c->i;
1368db02288aSLois Curfman McInnes     i_new[0] = -shift;
136902834360SBarry Smith     for (i=0; i<nrows; i++) {
1370a2744918SBarry Smith       ii    = starts[i];
1371a2744918SBarry Smith       lensi = lens[i];
1372a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1373a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
137402834360SBarry Smith       }
1375a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1376a2744918SBarry Smith       a_new      += lensi;
1377a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1378a2744918SBarry Smith       c->ilen[i]  = lensi;
137902834360SBarry Smith     }
13800452661fSBarry Smith     PetscFree(lens);
138102834360SBarry Smith   }
138202834360SBarry Smith   else {
138302834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
13840452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
138502834360SBarry Smith     ssmap = smap + shift;
138699141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1387cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
138817ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
138902834360SBarry Smith     /* determine lens of each row */
139002834360SBarry Smith     for (i=0; i<nrows; i++) {
1391d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
139202834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
139302834360SBarry Smith       lens[i] = 0;
139402834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1395d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
139602834360SBarry Smith           lens[i]++;
139702834360SBarry Smith         }
139802834360SBarry Smith       }
139902834360SBarry Smith     }
140017ab2063SBarry Smith     /* Create and fill new matrix */
1401a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
140299141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
140399141d43SSatish Balay 
1404e3372554SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(1,0,"");
140599141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1406e3372554SBarry Smith         SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros");
140799141d43SSatish Balay       }
140899141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
140908480c60SBarry Smith       C = *B;
141008480c60SBarry Smith     }
141108480c60SBarry Smith     else {
141202834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
141308480c60SBarry Smith     }
141499141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
141517ab2063SBarry Smith     for (i=0; i<nrows; i++) {
141699141d43SSatish Balay       row    = irow[i];
141717ab2063SBarry Smith       nznew  = 0;
141899141d43SSatish Balay       kstart = ai[row]+shift;
141999141d43SSatish Balay       kend   = kstart + a->ilen[row];
142099141d43SSatish Balay       mat_i  = c->i[i]+shift;
142199141d43SSatish Balay       mat_j  = c->j + mat_i;
142299141d43SSatish Balay       mat_a  = c->a + mat_i;
142399141d43SSatish Balay       mat_ilen = c->ilen + i;
142417ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
142599141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
142699141d43SSatish Balay           *mat_j++ = tcol - (!shift);
142799141d43SSatish Balay           *mat_a++ = a->a[k];
142899141d43SSatish Balay           (*mat_ilen)++;
142999141d43SSatish Balay 
143017ab2063SBarry Smith         }
143117ab2063SBarry Smith       }
143217ab2063SBarry Smith     }
143302834360SBarry Smith     /* Free work space */
143402834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
143599141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
143602834360SBarry Smith   }
14376d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
14386d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
143917ab2063SBarry Smith 
144017ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1441416022c9SBarry Smith   *B = C;
144217ab2063SBarry Smith   return 0;
144317ab2063SBarry Smith }
144417ab2063SBarry Smith 
1445a871dcd8SBarry Smith /*
144663b91edcSBarry Smith      note: This can only work for identity for row and col. It would
144763b91edcSBarry Smith    be good to check this and otherwise generate an error.
1448a871dcd8SBarry Smith */
14495615d1e5SSatish Balay #undef __FUNC__
14505615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
14518f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1452a871dcd8SBarry Smith {
145363b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
145408480c60SBarry Smith   int        ierr;
145563b91edcSBarry Smith   Mat        outA;
145663b91edcSBarry Smith 
1457e3372554SBarry Smith   if (fill != 0) SETERRQ(1,0,"Only fill=0 supported");
1458a871dcd8SBarry Smith 
145963b91edcSBarry Smith   outA          = inA;
146063b91edcSBarry Smith   inA->factor   = FACTOR_LU;
146163b91edcSBarry Smith   a->row        = row;
146263b91edcSBarry Smith   a->col        = col;
146363b91edcSBarry Smith 
146494a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
14650452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
146694a9d846SBarry Smith   }
146763b91edcSBarry Smith 
146808480c60SBarry Smith   if (!a->diag) {
146908480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
147063b91edcSBarry Smith   }
147163b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1472a871dcd8SBarry Smith   return 0;
1473a871dcd8SBarry Smith }
1474a871dcd8SBarry Smith 
1475f0b747eeSBarry Smith #include "pinclude/plapack.h"
14765615d1e5SSatish Balay #undef __FUNC__
14775615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
14788f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1479f0b747eeSBarry Smith {
1480f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1481f0b747eeSBarry Smith   int        one = 1;
1482f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1483f0b747eeSBarry Smith   PLogFlops(a->nz);
1484f0b747eeSBarry Smith   return 0;
1485f0b747eeSBarry Smith }
1486f0b747eeSBarry Smith 
14875615d1e5SSatish Balay #undef __FUNC__
14885615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
14898f6be9afSLois Curfman McInnes int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1490cddf8d76SBarry Smith                                     Mat **B)
1491cddf8d76SBarry Smith {
1492cddf8d76SBarry Smith   int ierr,i;
1493cddf8d76SBarry Smith 
1494cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
14950452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1496cddf8d76SBarry Smith   }
1497cddf8d76SBarry Smith 
1498cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1499905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1500cddf8d76SBarry Smith   }
1501cddf8d76SBarry Smith   return 0;
1502cddf8d76SBarry Smith }
1503cddf8d76SBarry Smith 
15045615d1e5SSatish Balay #undef __FUNC__
15055615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
15068f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
15075a838052SSatish Balay {
15085a838052SSatish Balay   *bs = 1;
15095a838052SSatish Balay   return 0;
15105a838052SSatish Balay }
15115a838052SSatish Balay 
15125615d1e5SSatish Balay #undef __FUNC__
15135615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15148f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15154dcbc457SBarry Smith {
1516e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
151706763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15188a047759SSatish Balay   int        start, end, *ai, *aj;
151906763907SSatish Balay   char       *table;
15208a047759SSatish Balay   shift = a->indexshift;
1521e4d965acSSatish Balay   m     = a->m;
1522e4d965acSSatish Balay   ai    = a->i;
15238a047759SSatish Balay   aj    = a->j+shift;
15248a047759SSatish Balay 
1525e3372554SBarry Smith   if (ov < 0)  SETERRQ(1,0,"illegal overlap value used");
152606763907SSatish Balay 
152706763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
152806763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
152906763907SSatish Balay 
1530e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1531b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1532e4d965acSSatish Balay     isz  = 0;
153306763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1534e4d965acSSatish Balay 
1535e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
15364dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
153777c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1538e4d965acSSatish Balay 
1539dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1540e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
154106763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
15424dcbc457SBarry Smith     }
154306763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
154406763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1545e4d965acSSatish Balay 
154604a348a9SBarry Smith     k = 0;
154704a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
154804a348a9SBarry Smith       n = isz;
154906763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1550e4d965acSSatish Balay         row   = nidx[k];
1551e4d965acSSatish Balay         start = ai[row];
1552e4d965acSSatish Balay         end   = ai[row+1];
155304a348a9SBarry Smith         for ( l = start; l<end ; l++){
15548a047759SSatish Balay           val = aj[l] + shift;
155506763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1556e4d965acSSatish Balay         }
1557e4d965acSSatish Balay       }
1558e4d965acSSatish Balay     }
1559029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1560e4d965acSSatish Balay   }
156104a348a9SBarry Smith   PetscFree(table);
156206763907SSatish Balay   PetscFree(nidx);
1563e4d965acSSatish Balay   return 0;
15644dcbc457SBarry Smith }
156517ab2063SBarry Smith 
15660513a670SBarry Smith /* -------------------------------------------------------------- */
15670513a670SBarry Smith #undef __FUNC__
15680513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
15690513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
15700513a670SBarry Smith {
15710513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
15720513a670SBarry Smith   Scalar     *vwork;
15730513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
15740513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
157556cd22aeSBarry Smith   IS         icolp,irowp;
15760513a670SBarry Smith 
157756cd22aeSBarry Smith   ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr);
157856cd22aeSBarry Smith   ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr);
157956cd22aeSBarry Smith   ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr);
158056cd22aeSBarry Smith   ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr);
15810513a670SBarry Smith 
15820513a670SBarry Smith   /* determine lengths of permuted rows */
15830513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
15840513a670SBarry Smith   for (i=0; i<m; i++ ) {
15850513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
15860513a670SBarry Smith   }
15870513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
15880513a670SBarry Smith   PetscFree(lens);
15890513a670SBarry Smith 
15900513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
15910513a670SBarry Smith   for (i=0; i<m; i++) {
15920513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15930513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
15940513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
15950513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15960513a670SBarry Smith   }
15970513a670SBarry Smith   PetscFree(cnew);
15980513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
15990513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
160056cd22aeSBarry Smith   ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr);
160156cd22aeSBarry Smith   ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr);
160256cd22aeSBarry Smith   ierr = ISDestroy(irowp); CHKERRQ(ierr);
160356cd22aeSBarry Smith   ierr = ISDestroy(icolp); CHKERRQ(ierr);
16040513a670SBarry Smith   return 0;
16050513a670SBarry Smith }
16060513a670SBarry Smith 
16075615d1e5SSatish Balay #undef __FUNC__
16085615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1609682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1610682d7d0cSBarry Smith {
1611682d7d0cSBarry Smith   static int called = 0;
1612682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1613682d7d0cSBarry Smith 
1614682d7d0cSBarry Smith   if (called) return 0; else called = 1;
161577c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
16160f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
16170f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
16180f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_no_inode: Do not use inodes\n");
16190f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1620682d7d0cSBarry Smith #if defined(HAVE_ESSL)
16210f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1622682d7d0cSBarry Smith #endif
1623682d7d0cSBarry Smith   return 0;
1624682d7d0cSBarry Smith }
16258f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1626a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1627a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1628a93ec695SBarry Smith 
1629682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
163017ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
163117ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1632416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1633416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
163417ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
163517ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
163617ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
163717ab2063SBarry Smith        MatRelax_SeqAIJ,
163817ab2063SBarry Smith        MatTranspose_SeqAIJ,
16397264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1640f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
164117ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
164217ab2063SBarry Smith        MatCompress_SeqAIJ,
164317ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
164417ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
164517ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
164617ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
164794a9d846SBarry Smith        0,0,
16483d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1649cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
16507eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1651682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1652f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
16535a838052SSatish Balay        MatScale_SeqAIJ,0,0,
16546945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
16556945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
16563b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
16573b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
16583b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1659a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1660a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
16610513a670SBarry Smith        MatColoringPatch_SeqAIJ,
16620513a670SBarry Smith        0,
16630513a670SBarry Smith        MatPermute_SeqAIJ};
166417ab2063SBarry Smith 
166517ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
166617ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
166717ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
166817ab2063SBarry Smith 
16695615d1e5SSatish Balay #undef __FUNC__
16705615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
167117ab2063SBarry Smith /*@C
1672682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
16730d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
16746e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
16752bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
16762bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
167717ab2063SBarry Smith 
167817ab2063SBarry Smith    Input Parameters:
1679029af93fSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
168017ab2063SBarry Smith .  m - number of rows
168117ab2063SBarry Smith .  n - number of columns
168217ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
16832bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
16842bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
168517ab2063SBarry Smith 
168617ab2063SBarry Smith    Output Parameter:
1687416022c9SBarry Smith .  A - the matrix
168817ab2063SBarry Smith 
168917ab2063SBarry Smith    Notes:
169017ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
169117ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
16920002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
169344cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
169417ab2063SBarry Smith 
169517ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1696a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
16973d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
16986da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
169917ab2063SBarry Smith 
1700682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
17014fca80b9SLois Curfman McInnes    improve numerical efficiency of matrix-vector products and solves. We
1702682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
17036c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
17046c7ebb05SLois Curfman McInnes 
17056c7ebb05SLois Curfman McInnes    Options Database Keys:
17066c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
17076e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
17086e62573dSLois Curfman McInnes $        (max limit=5)
17096e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
17106e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
17116e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
171217ab2063SBarry Smith 
171317ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
171417ab2063SBarry Smith @*/
1715416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
171617ab2063SBarry Smith {
1717416022c9SBarry Smith   Mat        B;
1718416022c9SBarry Smith   Mat_SeqAIJ *b;
17196945ee14SBarry Smith   int        i, len, ierr, flg,size;
17206945ee14SBarry Smith 
17216945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1722e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
1723d5d45c9bSBarry Smith 
1724416022c9SBarry Smith   *A                  = 0;
1725d4bb536fSBarry Smith   PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView);
1726416022c9SBarry Smith   PLogObjectCreate(B);
17270452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
172844cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1729416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1730416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1731416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1732416022c9SBarry Smith   B->factor           = 0;
1733416022c9SBarry Smith   B->lupivotthreshold = 1.0;
173490f02eecSBarry Smith   B->mapping          = 0;
17357a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
173669957df2SSatish Balay                           &flg); CHKERRQ(ierr);
17377a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
17387a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
17397a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1740416022c9SBarry Smith   b->row              = 0;
1741416022c9SBarry Smith   b->col              = 0;
1742416022c9SBarry Smith   b->indexshift       = 0;
1743b810aeb4SBarry Smith   b->reallocs         = 0;
174469957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
174569957df2SSatish Balay   if (flg) b->indexshift = -1;
174617ab2063SBarry Smith 
174744cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
174844cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
17490452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1750b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
17517b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
17527b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1753416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
175417ab2063SBarry Smith     nz = nz*m;
175517ab2063SBarry Smith   }
175617ab2063SBarry Smith   else {
175717ab2063SBarry Smith     nz = 0;
1758416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
175917ab2063SBarry Smith   }
176017ab2063SBarry Smith 
176117ab2063SBarry Smith   /* allocate the matrix space */
1762416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
17630452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1764416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1765cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1766416022c9SBarry Smith   b->i  = b->j + nz;
1767416022c9SBarry Smith   b->singlemalloc = 1;
176817ab2063SBarry Smith 
1769416022c9SBarry Smith   b->i[0] = -b->indexshift;
177017ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1771416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
177217ab2063SBarry Smith   }
177317ab2063SBarry Smith 
1774416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
17750452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1776f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1777416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
177817ab2063SBarry Smith 
1779416022c9SBarry Smith   b->nz               = 0;
1780416022c9SBarry Smith   b->maxnz            = nz;
1781416022c9SBarry Smith   b->sorted           = 0;
1782416022c9SBarry Smith   b->roworiented      = 1;
1783416022c9SBarry Smith   b->nonew            = 0;
1784416022c9SBarry Smith   b->diag             = 0;
1785416022c9SBarry Smith   b->solve_work       = 0;
178671bd300dSLois Curfman McInnes   b->spptr            = 0;
1787754ec7b1SSatish Balay   b->inode.node_count = 0;
1788754ec7b1SSatish Balay   b->inode.size       = 0;
17896c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
17906c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
17914e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
179217ab2063SBarry Smith 
1793416022c9SBarry Smith   *A = B;
17944e220ebcSLois Curfman McInnes 
17954b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
17964b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
179769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
179869957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
17994b14c69eSBarry Smith #endif
180069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
180169957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
180269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
180369957df2SSatish Balay   if (flg) {
1804e3372554SBarry Smith     if (!b->indexshift) SETERRQ(1,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1805416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
180617ab2063SBarry Smith   }
180769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
180869957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
180917ab2063SBarry Smith   return 0;
181017ab2063SBarry Smith }
181117ab2063SBarry Smith 
18125615d1e5SSatish Balay #undef __FUNC__
18135615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
18143d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
181517ab2063SBarry Smith {
1816416022c9SBarry Smith   Mat        C;
1817416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
181808480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
181917ab2063SBarry Smith 
18204043dd9cSLois Curfman McInnes   *B = 0;
1821d4bb536fSBarry Smith   PetscHeaderCreate(C,_p_Mat,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView);
1822416022c9SBarry Smith   PLogObjectCreate(C);
18230452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
182441c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1825416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1826416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1827416022c9SBarry Smith   C->factor     = A->factor;
1828416022c9SBarry Smith   c->row        = 0;
1829416022c9SBarry Smith   c->col        = 0;
1830416022c9SBarry Smith   c->indexshift = shift;
1831c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
183217ab2063SBarry Smith 
183344cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
183444cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
183544cd7ae7SLois Curfman McInnes   C->M          = a->m;
183644cd7ae7SLois Curfman McInnes   C->N          = a->n;
183717ab2063SBarry Smith 
18380452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
18390452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
184017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1841416022c9SBarry Smith     c->imax[i] = a->imax[i];
1842416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
184317ab2063SBarry Smith   }
184417ab2063SBarry Smith 
184517ab2063SBarry Smith   /* allocate the matrix space */
1846416022c9SBarry Smith   c->singlemalloc = 1;
1847416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
18480452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1849416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1850416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1851416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
185217ab2063SBarry Smith   if (m > 0) {
1853416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
185408480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1855416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
185617ab2063SBarry Smith     }
185708480c60SBarry Smith   }
185817ab2063SBarry Smith 
1859f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1860416022c9SBarry Smith   c->sorted      = a->sorted;
1861416022c9SBarry Smith   c->roworiented = a->roworiented;
1862416022c9SBarry Smith   c->nonew       = a->nonew;
18637a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
186417ab2063SBarry Smith 
1865416022c9SBarry Smith   if (a->diag) {
18660452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1867416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
186817ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1869416022c9SBarry Smith       c->diag[i] = a->diag[i];
187017ab2063SBarry Smith     }
187117ab2063SBarry Smith   }
1872416022c9SBarry Smith   else c->diag          = 0;
18736c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
18746c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1875754ec7b1SSatish Balay   if (a->inode.size){
1876daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1877754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1878daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1879754ec7b1SSatish Balay   } else {
1880754ec7b1SSatish Balay     c->inode.size       = 0;
1881754ec7b1SSatish Balay     c->inode.node_count = 0;
1882754ec7b1SSatish Balay   }
1883416022c9SBarry Smith   c->nz                 = a->nz;
1884416022c9SBarry Smith   c->maxnz              = a->maxnz;
1885416022c9SBarry Smith   c->solve_work         = 0;
188676dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1887754ec7b1SSatish Balay 
1888416022c9SBarry Smith   *B = C;
188917ab2063SBarry Smith   return 0;
189017ab2063SBarry Smith }
189117ab2063SBarry Smith 
18925615d1e5SSatish Balay #undef __FUNC__
18935615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
189419bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
189517ab2063SBarry Smith {
1896416022c9SBarry Smith   Mat_SeqAIJ   *a;
1897416022c9SBarry Smith   Mat          B;
189817699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1899bcd2baecSBarry Smith   MPI_Comm     comm;
190017ab2063SBarry Smith 
190119bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
190217699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
1903e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"view must have one processor");
190490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
190577c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1906e3372554SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object in file");
190717ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
190817ab2063SBarry Smith 
190917ab2063SBarry Smith   /* read in row lengths */
19100452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
191177c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
191217ab2063SBarry Smith 
191317ab2063SBarry Smith   /* create our matrix */
1914416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1915416022c9SBarry Smith   B = *A;
1916416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1917416022c9SBarry Smith   shift = a->indexshift;
191817ab2063SBarry Smith 
191917ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
192077c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
192117ab2063SBarry Smith   if (shift) {
192217ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1923416022c9SBarry Smith       a->j[i] += 1;
192417ab2063SBarry Smith     }
192517ab2063SBarry Smith   }
192617ab2063SBarry Smith 
192717ab2063SBarry Smith   /* read in nonzero values */
192877c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
192917ab2063SBarry Smith 
193017ab2063SBarry Smith   /* set matrix "i" values */
1931416022c9SBarry Smith   a->i[0] = -shift;
193217ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1933416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1934416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
193517ab2063SBarry Smith   }
19360452661fSBarry Smith   PetscFree(rowlengths);
193717ab2063SBarry Smith 
19386d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
19396d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
194017ab2063SBarry Smith   return 0;
194117ab2063SBarry Smith }
194217ab2063SBarry Smith 
19435615d1e5SSatish Balay #undef __FUNC__
19445615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
19458f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
19467264ac53SSatish Balay {
19477264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
19487264ac53SSatish Balay 
1949e3372554SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,0,"Matrices must be same type");
19507264ac53SSatish Balay 
19517264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
19527264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1953bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
195477c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1955bcd2baecSBarry Smith   }
19567264ac53SSatish Balay 
19577264ac53SSatish Balay   /* if the a->i are the same */
19588108c231SLois Curfman McInnes   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
195977c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
19607264ac53SSatish Balay   }
19617264ac53SSatish Balay 
19627264ac53SSatish Balay   /* if a->j are the same */
1963bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
196477c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1965bcd2baecSBarry Smith   }
1966bcd2baecSBarry Smith 
1967bcd2baecSBarry Smith   /* if a->a are the same */
196819bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
196977c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
19707264ac53SSatish Balay   }
197177c4ece6SBarry Smith   *flg = PETSC_TRUE;
19727264ac53SSatish Balay   return 0;
19737264ac53SSatish Balay 
19747264ac53SSatish Balay }
1975