xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 2e44a96cb8a880801236bd5ac15638be8ae545a8)
117ab2063SBarry Smith #ifndef lint
2*2e44a96cSLois Curfman McInnes static char vcid[] = "$Id: aij.c,v 1.224 1997/06/19 00:02:06 curfman Exp curfman $";
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"
15e4d965acSSatish Balay #include "petsc.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) {
336416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3374e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
3384e220ebcSLois Curfman McInnes     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz);
33917ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
34017ab2063SBarry Smith 
34117ab2063SBarry Smith     for (i=0; i<m; i++) {
342416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
34317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
3446945ee14SBarry Smith         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]),
345416022c9SBarry Smith                    imag(a->a[j]));
34617ab2063SBarry Smith #else
3477a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
34817ab2063SBarry Smith #endif
34917ab2063SBarry Smith       }
35017ab2063SBarry Smith     }
35117ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
35217ab2063SBarry Smith   }
353a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
35444cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
35544cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
35644cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
35744cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
358766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0 && real(a->a[j]) != 0.0)
35944cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
360766eeae4SLois Curfman McInnes         else if (imag(a->a[j]) < 0.0 && real(a->a[j]) != 0.0)
361766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
36244cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
36344cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
36444cd7ae7SLois Curfman McInnes #else
36544cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
36644cd7ae7SLois Curfman McInnes #endif
36744cd7ae7SLois Curfman McInnes       }
36844cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
36944cd7ae7SLois Curfman McInnes     }
37044cd7ae7SLois Curfman McInnes   }
371496be53dSLois Curfman McInnes   else if (format == VIEWER_FORMAT_ASCII_SYMMODU) {
372496be53dSLois Curfman McInnes     int nzd=0, fshift=1, *sptr;
373*2e44a96cSLois Curfman McInnes     sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr);
374496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
375496be53dSLois Curfman McInnes       sptr[i] = nzd+1;
376496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
377496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
378496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX)
379496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) nzd++;
380496be53dSLois Curfman McInnes #else
381496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) nzd++;
382496be53dSLois Curfman McInnes #endif
383496be53dSLois Curfman McInnes         }
384496be53dSLois Curfman McInnes       }
385496be53dSLois Curfman McInnes     }
386*2e44a96cSLois Curfman McInnes     sptr[m] = nzd+1;
387496be53dSLois Curfman McInnes     fprintf(fd," %d %d\n\n",m,nzd);
388*2e44a96cSLois Curfman McInnes     for ( i=0; i<m+1; i+=6 ) {
389*2e44a96cSLois 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]);
390*2e44a96cSLois 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]);
391*2e44a96cSLois 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]);
392*2e44a96cSLois Curfman McInnes       else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]);
393*2e44a96cSLois Curfman McInnes       else if (i<m)   fprintf(fd," %d %d\n",sptr[i],sptr[i+1]);
3947272d637SLois Curfman McInnes       else            fprintf(fd," %d\n",sptr[i]);
395496be53dSLois Curfman McInnes     }
396496be53dSLois Curfman McInnes     fprintf(fd,"\n");
397496be53dSLois Curfman McInnes     PetscFree(sptr);
398496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
399496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
400496be53dSLois Curfman McInnes         if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift);
401496be53dSLois Curfman McInnes       }
402496be53dSLois Curfman McInnes       fprintf(fd,"\n");
403496be53dSLois Curfman McInnes     }
404496be53dSLois Curfman McInnes     fprintf(fd,"\n");
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) {
408496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX)
409496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0)
410496be53dSLois Curfman McInnes             fprintf(fd," %18.16e %18.16e ",real(a->a[j]),imag(a->a[j]));
411496be53dSLois Curfman McInnes #else
412496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]);
413496be53dSLois Curfman McInnes #endif
414496be53dSLois Curfman McInnes         }
415496be53dSLois Curfman McInnes       }
416496be53dSLois Curfman McInnes       fprintf(fd,"\n");
417496be53dSLois Curfman McInnes     }
418496be53dSLois Curfman McInnes   }
41917ab2063SBarry Smith   else {
42017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
42117ab2063SBarry Smith       fprintf(fd,"row %d:",i);
422416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
42317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
424766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0) {
425416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
426766eeae4SLois Curfman McInnes         } else if (imag(a->a[j]) < 0.0) {
427766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
42817ab2063SBarry Smith         }
42917ab2063SBarry Smith         else {
430416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
43117ab2063SBarry Smith         }
43217ab2063SBarry Smith #else
433416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
43417ab2063SBarry Smith #endif
43517ab2063SBarry Smith       }
43617ab2063SBarry Smith       fprintf(fd,"\n");
43717ab2063SBarry Smith     }
43817ab2063SBarry Smith   }
43917ab2063SBarry Smith   fflush(fd);
440416022c9SBarry Smith   return 0;
441416022c9SBarry Smith }
442416022c9SBarry Smith 
4435615d1e5SSatish Balay #undef __FUNC__
4445615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Draw"
4458f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
446416022c9SBarry Smith {
447416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
448cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
4490513a670SBarry Smith   int         format;
45094a9d846SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r,maxv = 0.0;
451bcd2baecSBarry Smith   Draw        draw;
452cddf8d76SBarry Smith   DrawButton  button;
45319bcc07fSBarry Smith   PetscTruth  isnull;
454cddf8d76SBarry Smith 
4550513a670SBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
4560513a670SBarry Smith   ierr = DrawClear(draw); CHKERRQ(ierr);
4570513a670SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
45819bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
45919bcc07fSBarry Smith 
460416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
461416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
462416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
463416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
4640513a670SBarry Smith 
4650513a670SBarry Smith   if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4660513a670SBarry Smith     /* Blue for negative, Cyan for zero and  Red for positive */
467cddf8d76SBarry Smith     color = DRAW_BLUE;
468416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
469cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
470416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
471cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
472cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
473cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
474cddf8d76SBarry Smith #else
475cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
476cddf8d76SBarry Smith #endif
477cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
478cddf8d76SBarry Smith       }
479cddf8d76SBarry Smith     }
480cddf8d76SBarry Smith     color = DRAW_CYAN;
481cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
482cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
483cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
484cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
485cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
486cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
487cddf8d76SBarry Smith       }
488cddf8d76SBarry Smith     }
489cddf8d76SBarry Smith     color = DRAW_RED;
490cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
491cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
492cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
493cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
494cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
495cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
496cddf8d76SBarry Smith #else
497cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
498cddf8d76SBarry Smith #endif
499cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
500416022c9SBarry Smith       }
501416022c9SBarry Smith     }
5020513a670SBarry Smith   } else {
5030513a670SBarry Smith     /* use contour shading to indicate magnitude of values */
5040513a670SBarry Smith     /* first determine max of all nonzero values */
5050513a670SBarry Smith     int    nz = a->nz,count;
5060513a670SBarry Smith     Draw   popup;
5070513a670SBarry Smith 
5080513a670SBarry Smith     for ( i=0; i<nz; i++ ) {
5090513a670SBarry Smith       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
5100513a670SBarry Smith     }
5110513a670SBarry Smith     ierr = DrawCreatePopUp(draw,&popup); CHKERRQ(ierr);
5120513a670SBarry Smith     ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr);
5130513a670SBarry Smith     count = 0;
5140513a670SBarry Smith     for ( i=0; i<m; i++ ) {
5150513a670SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
5160513a670SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5170513a670SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5180513a670SBarry Smith         color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5190513a670SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5200513a670SBarry Smith         count++;
5210513a670SBarry Smith       }
5220513a670SBarry Smith     }
5230513a670SBarry Smith   }
524416022c9SBarry Smith   DrawFlush(draw);
525cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
526cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
527cddf8d76SBarry Smith 
528cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
5296945ee14SBarry Smith   ierr = DrawCheckResizedWindow(draw);
530cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
531cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
532cddf8d76SBarry Smith     DrawClear(draw);
533cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
534cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
535cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
536cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
537cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
538cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
539cddf8d76SBarry Smith     w *= scale; h *= scale;
540cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5410513a670SBarry Smith     if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
5420513a670SBarry Smith       /* Blue for negative, Cyan for zero and  Red for positive */
543cddf8d76SBarry Smith       color = DRAW_BLUE;
544cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
545cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
546cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
547cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
548cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
549cddf8d76SBarry Smith           if (real(a->a[j]) >=  0.) continue;
550cddf8d76SBarry Smith #else
551cddf8d76SBarry Smith           if (a->a[j] >=  0.) continue;
552cddf8d76SBarry Smith #endif
553cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
554cddf8d76SBarry Smith         }
555cddf8d76SBarry Smith       }
556cddf8d76SBarry Smith       color = DRAW_CYAN;
557cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
558cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
559cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
560cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
561cddf8d76SBarry Smith           if (a->a[j] !=  0.) continue;
562cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
563cddf8d76SBarry Smith         }
564cddf8d76SBarry Smith       }
565cddf8d76SBarry Smith       color = DRAW_RED;
566cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
567cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
568cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
569cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
570cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
571cddf8d76SBarry Smith           if (real(a->a[j]) <=  0.) continue;
572cddf8d76SBarry Smith #else
573cddf8d76SBarry Smith           if (a->a[j] <=  0.) continue;
574cddf8d76SBarry Smith #endif
575cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
576cddf8d76SBarry Smith         }
577cddf8d76SBarry Smith       }
5780513a670SBarry Smith     } else {
5790513a670SBarry Smith       /* use contour shading to indicate magnitude of values */
5800513a670SBarry Smith       int count = 0;
5810513a670SBarry Smith       for ( i=0; i<m; i++ ) {
5820513a670SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
5830513a670SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5840513a670SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
5850513a670SBarry Smith           color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5860513a670SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5870513a670SBarry Smith           count++;
5880513a670SBarry Smith         }
5890513a670SBarry Smith       }
5900513a670SBarry Smith     }
5910513a670SBarry Smith 
5926945ee14SBarry Smith     ierr = DrawCheckResizedWindow(draw);
593cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
594cddf8d76SBarry Smith   }
595416022c9SBarry Smith   return 0;
596416022c9SBarry Smith }
597416022c9SBarry Smith 
5985615d1e5SSatish Balay #undef __FUNC__
599c22c1629SBarry Smith #define __FUNC__ "MatView_SeqAIJ" /* ADIC Ignore */
6008f6be9afSLois Curfman McInnes int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
601416022c9SBarry Smith {
602416022c9SBarry Smith   Mat         A = (Mat) obj;
603416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
604bcd2baecSBarry Smith   ViewerType  vtype;
605bcd2baecSBarry Smith   int         ierr;
606416022c9SBarry Smith 
607bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
608bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
609416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
610416022c9SBarry Smith   }
611bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
612416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
613416022c9SBarry Smith   }
614bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
615416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
616416022c9SBarry Smith   }
617bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
618bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
61917ab2063SBarry Smith   }
62017ab2063SBarry Smith   return 0;
62117ab2063SBarry Smith }
62219bcc07fSBarry Smith 
623c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
6245615d1e5SSatish Balay #undef __FUNC__
6255615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ"
6268f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
62717ab2063SBarry Smith {
628416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
62941c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
63094a9d846SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax;
631416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
63217ab2063SBarry Smith 
6336d4a8577SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
63417ab2063SBarry Smith 
63594a9d846SBarry Smith   rmax = ailen[0]; /* determine row with most nonzeros */
63617ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
637416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
63817ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
63994a9d846SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
64017ab2063SBarry Smith     if (fshift) {
641416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
64217ab2063SBarry Smith       N = ailen[i];
64317ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
64417ab2063SBarry Smith         ip[j-fshift] = ip[j];
64517ab2063SBarry Smith         ap[j-fshift] = ap[j];
64617ab2063SBarry Smith       }
64717ab2063SBarry Smith     }
64817ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
64917ab2063SBarry Smith   }
65017ab2063SBarry Smith   if (m) {
65117ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
65217ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
65317ab2063SBarry Smith   }
65417ab2063SBarry Smith   /* reset ilen and imax for each row */
65517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
65617ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
65717ab2063SBarry Smith   }
658416022c9SBarry Smith   a->nz = ai[m] + shift;
65917ab2063SBarry Smith 
66017ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
661416022c9SBarry Smith   if (fshift && a->diag) {
6620452661fSBarry Smith     PetscFree(a->diag);
663416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
664416022c9SBarry Smith     a->diag = 0;
66517ab2063SBarry Smith   }
6664e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
6674e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
6684e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
669b810aeb4SBarry Smith            a->reallocs);
67094a9d846SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax);
671dd5f02e7SSatish Balay   a->reallocs          = 0;
6724e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
6734e220ebcSLois Curfman McInnes 
67476dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
67541c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
67617ab2063SBarry Smith   return 0;
67717ab2063SBarry Smith }
67817ab2063SBarry Smith 
6795615d1e5SSatish Balay #undef __FUNC__
6805615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ"
6818f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A)
68217ab2063SBarry Smith {
683416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
684cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
68517ab2063SBarry Smith   return 0;
68617ab2063SBarry Smith }
687416022c9SBarry Smith 
6885615d1e5SSatish Balay #undef __FUNC__
6895615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ"
69017ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
69117ab2063SBarry Smith {
692416022c9SBarry Smith   Mat        A  = (Mat) obj;
693416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
694d5d45c9bSBarry Smith 
69517ab2063SBarry Smith #if defined(PETSC_LOG)
696416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
69717ab2063SBarry Smith #endif
6980452661fSBarry Smith   PetscFree(a->a);
6990452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
7000452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
7010452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
7020452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
7030452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
70476dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
7050452661fSBarry Smith   PetscFree(a);
706eed86810SBarry Smith 
707f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
708f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
70917ab2063SBarry Smith   return 0;
71017ab2063SBarry Smith }
71117ab2063SBarry Smith 
7125615d1e5SSatish Balay #undef __FUNC__
7135615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
7148f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A)
71517ab2063SBarry Smith {
71617ab2063SBarry Smith   return 0;
71717ab2063SBarry Smith }
71817ab2063SBarry Smith 
7195615d1e5SSatish Balay #undef __FUNC__
7205615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
7218f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op)
72217ab2063SBarry Smith {
723416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7246d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
7256d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
7266d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
727219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
7286d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
729c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
73096854ed6SLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
7316d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
7326d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
733219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7346d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7356d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
73690f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
7372b362799SSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES)
73894a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
7396d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
740e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");}
7416d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
7426d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
7436d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
7446d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
7456d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
746e2f28af5SBarry Smith   else
747e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
74817ab2063SBarry Smith   return 0;
74917ab2063SBarry Smith }
75017ab2063SBarry Smith 
7515615d1e5SSatish Balay #undef __FUNC__
7525615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
7538f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
75417ab2063SBarry Smith {
755416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
756416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
75717ab2063SBarry Smith   Scalar     *x, zero = 0.0;
75817ab2063SBarry Smith 
75917ab2063SBarry Smith   VecSet(&zero,v);
76090f02eecSBarry Smith   VecGetArray_Fast(v,x); VecGetLocalSize(v,&n);
761e3372554SBarry Smith   if (n != a->m) SETERRQ(1,0,"Nonconforming matrix and vector");
762416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
763416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
764416022c9SBarry Smith       if (a->j[j]+shift == i) {
765416022c9SBarry Smith         x[i] = a->a[j];
76617ab2063SBarry Smith         break;
76717ab2063SBarry Smith       }
76817ab2063SBarry Smith     }
76917ab2063SBarry Smith   }
77017ab2063SBarry Smith   return 0;
77117ab2063SBarry Smith }
77217ab2063SBarry Smith 
77317ab2063SBarry Smith /* -------------------------------------------------------*/
77417ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
77517ab2063SBarry Smith /* -------------------------------------------------------*/
7765615d1e5SSatish Balay #undef __FUNC__
7775615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
77844cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
77917ab2063SBarry Smith {
780416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
78117ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
782416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
78317ab2063SBarry Smith 
78490f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
785cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
78617ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
78717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
788416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
789416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
790416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
79117ab2063SBarry Smith     alpha = x[i];
79217ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
79317ab2063SBarry Smith   }
794416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
79517ab2063SBarry Smith   return 0;
79617ab2063SBarry Smith }
797d5d45c9bSBarry Smith 
7985615d1e5SSatish Balay #undef __FUNC__
7995615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
80044cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
80117ab2063SBarry Smith {
802416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
80317ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
804416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
80517ab2063SBarry Smith 
80690f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
80717ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
80817ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
80917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
810416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
811416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
812416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
81317ab2063SBarry Smith     alpha = x[i];
81417ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
81517ab2063SBarry Smith   }
81690f02eecSBarry Smith   PLogFlops(2*a->nz);
81717ab2063SBarry Smith   return 0;
81817ab2063SBarry Smith }
81917ab2063SBarry Smith 
8205615d1e5SSatish Balay #undef __FUNC__
8215615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
82244cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
82317ab2063SBarry Smith {
824416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
82517ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
8269ea0dfa2SSatish Balay   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j;
82717ab2063SBarry Smith 
82890f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
82917ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
830416022c9SBarry Smith   idx  = a->j;
831416022c9SBarry Smith   v    = a->a;
832416022c9SBarry Smith   ii   = a->i;
83317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8349ea0dfa2SSatish Balay     jrow = ii[i];
8359ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
83617ab2063SBarry Smith     sum  = 0.0;
8379ea0dfa2SSatish Balay     /* while (n--) sum += *v++ * x[*idx++]; */
8389ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8399ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8409ea0dfa2SSatish Balay      }
84117ab2063SBarry Smith     y[i] = sum;
84217ab2063SBarry Smith   }
843416022c9SBarry Smith   PLogFlops(2*a->nz - m);
84417ab2063SBarry Smith   return 0;
84517ab2063SBarry Smith }
84617ab2063SBarry Smith 
8475615d1e5SSatish Balay #undef __FUNC__
8485615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
84944cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
85017ab2063SBarry Smith {
851416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
85217ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
8539ea0dfa2SSatish Balay   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j;
8549ea0dfa2SSatish Balay 
85517ab2063SBarry Smith 
85690f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); VecGetArray_Fast(zz,z);
85717ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
858cddf8d76SBarry Smith   idx  = a->j;
859cddf8d76SBarry Smith   v    = a->a;
860cddf8d76SBarry Smith   ii   = a->i;
86117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8629ea0dfa2SSatish Balay     jrow = ii[i];
8639ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
86417ab2063SBarry Smith     sum  = y[i];
8659ea0dfa2SSatish Balay     /* while (n--) sum += *v++ * x[*idx++]; */
8669ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8679ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8689ea0dfa2SSatish Balay      }
86917ab2063SBarry Smith     z[i] = sum;
87017ab2063SBarry Smith   }
871416022c9SBarry Smith   PLogFlops(2*a->nz);
87217ab2063SBarry Smith   return 0;
87317ab2063SBarry Smith }
87417ab2063SBarry Smith 
87517ab2063SBarry Smith /*
87617ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
87717ab2063SBarry Smith */
87817ab2063SBarry Smith 
8795615d1e5SSatish Balay #undef __FUNC__
8805615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
881416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
88217ab2063SBarry Smith {
883416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
884416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
88517ab2063SBarry Smith 
8860452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
887416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
888416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
889416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
890416022c9SBarry Smith       if (a->j[j]+shift == i) {
89117ab2063SBarry Smith         diag[i] = j - shift;
89217ab2063SBarry Smith         break;
89317ab2063SBarry Smith       }
89417ab2063SBarry Smith     }
89517ab2063SBarry Smith   }
896416022c9SBarry Smith   a->diag = diag;
89717ab2063SBarry Smith   return 0;
89817ab2063SBarry Smith }
89917ab2063SBarry Smith 
9005615d1e5SSatish Balay #undef __FUNC__
9015615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
90244cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
90317ab2063SBarry Smith                            double fshift,int its,Vec xx)
90417ab2063SBarry Smith {
905416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
906416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
907d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
90817ab2063SBarry Smith 
90990f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(bb,b);
910416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
911416022c9SBarry Smith   diag = a->diag;
912416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
91317ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
91417ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
91517ab2063SBarry Smith     bs = b + shift;
91617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
917416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
918416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
919416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
920416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
92117ab2063SBarry Smith         sum  = b[i]*d/omega;
92217ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
92317ab2063SBarry Smith         x[i] = sum;
92417ab2063SBarry Smith     }
92517ab2063SBarry Smith     return 0;
92617ab2063SBarry Smith   }
92717ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
928e3372554SBarry Smith     SETERRQ(1,0,"SOR_APPLY_LOWER is not done");
92917ab2063SBarry Smith   }
930416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
93117ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
93217ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
93317ab2063SBarry Smith 
93417ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
93517ab2063SBarry Smith 
93617ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
93717ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
93817ab2063SBarry Smith     is the relaxation factor.
93917ab2063SBarry Smith     */
9400452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
94117ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
94217ab2063SBarry Smith 
94317ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
94417ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
945416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
946416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
947416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
948416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
94917ab2063SBarry Smith       sum  = b[i];
95017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
95117ab2063SBarry Smith       x[i] = omega*(sum/d);
95217ab2063SBarry Smith     }
95317ab2063SBarry Smith 
95417ab2063SBarry Smith     /*  t = b - (2*E - D)x */
955416022c9SBarry Smith     v = a->a;
95617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
95717ab2063SBarry Smith 
95817ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
959416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
960416022c9SBarry Smith     diag = a->diag;
96117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
962416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
963416022c9SBarry Smith       n    = diag[i] - a->i[i];
964416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
965416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
96617ab2063SBarry Smith       sum  = t[i];
96717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
96817ab2063SBarry Smith       t[i] = omega*(sum/d);
96917ab2063SBarry Smith     }
97017ab2063SBarry Smith 
97117ab2063SBarry Smith     /*  x = x + t */
97217ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
9730452661fSBarry Smith     PetscFree(t);
97417ab2063SBarry Smith     return 0;
97517ab2063SBarry Smith   }
97617ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
97717ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
97817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
979416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
980416022c9SBarry Smith         n    = diag[i] - a->i[i];
981416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
982416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
98317ab2063SBarry Smith         sum  = b[i];
98417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
98517ab2063SBarry Smith         x[i] = omega*(sum/d);
98617ab2063SBarry Smith       }
98717ab2063SBarry Smith       xb = x;
98817ab2063SBarry Smith     }
98917ab2063SBarry Smith     else xb = b;
99017ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
99117ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
99217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
993416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
99417ab2063SBarry Smith       }
99517ab2063SBarry Smith     }
99617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
99717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
998416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
999416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
1000416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1001416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
100217ab2063SBarry Smith         sum  = xb[i];
100317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
100417ab2063SBarry Smith         x[i] = omega*(sum/d);
100517ab2063SBarry Smith       }
100617ab2063SBarry Smith     }
100717ab2063SBarry Smith     its--;
100817ab2063SBarry Smith   }
100917ab2063SBarry Smith   while (its--) {
101017ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
101117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1012416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1013416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1014416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1015416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
101617ab2063SBarry Smith         sum  = b[i];
101717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10187e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
101917ab2063SBarry Smith       }
102017ab2063SBarry Smith     }
102117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
102217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1023416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1024416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1025416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1026416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
102717ab2063SBarry Smith         sum  = b[i];
102817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10297e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
103017ab2063SBarry Smith       }
103117ab2063SBarry Smith     }
103217ab2063SBarry Smith   }
103317ab2063SBarry Smith   return 0;
103417ab2063SBarry Smith }
103517ab2063SBarry Smith 
10365615d1e5SSatish Balay #undef __FUNC__
10375615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10388f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
103917ab2063SBarry Smith {
1040416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
10414e220ebcSLois Curfman McInnes 
10424e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
10434e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
10444e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
10454e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
10464e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10474e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
10484e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
10494e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
10504e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
10514e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
10524e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
10534e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
10544e220ebcSLois Curfman McInnes   info->memory         = A->mem;
10554e220ebcSLois Curfman McInnes   if (A->factor) {
10564e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10574e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
10584e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
10594e220ebcSLois Curfman McInnes   } else {
10604e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
10614e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
10624e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
10634e220ebcSLois Curfman McInnes   }
106417ab2063SBarry Smith   return 0;
106517ab2063SBarry Smith }
106617ab2063SBarry Smith 
106717ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
106817ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
106917ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
107017ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
107117ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
107217ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
107317ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
107417ab2063SBarry Smith 
10755615d1e5SSatish Balay #undef __FUNC__
10765615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
10778f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
107817ab2063SBarry Smith {
1079416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1080416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
108117ab2063SBarry Smith 
108277c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
108317ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
108417ab2063SBarry Smith   if (diag) {
108517ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1086e3372554SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range");
1087416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1088416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1089416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1090416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
109117ab2063SBarry Smith       }
109217ab2063SBarry Smith       else {
109317ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
109417ab2063SBarry Smith         CHKERRQ(ierr);
109517ab2063SBarry Smith       }
109617ab2063SBarry Smith     }
109717ab2063SBarry Smith   }
109817ab2063SBarry Smith   else {
109917ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1100e3372554SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range");
1101416022c9SBarry Smith       a->ilen[rows[i]] = 0;
110217ab2063SBarry Smith     }
110317ab2063SBarry Smith   }
110417ab2063SBarry Smith   ISRestoreIndices(is,&rows);
110543a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
110617ab2063SBarry Smith   return 0;
110717ab2063SBarry Smith }
110817ab2063SBarry Smith 
11095615d1e5SSatish Balay #undef __FUNC__
11105615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11118f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
111217ab2063SBarry Smith {
1113416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1114416022c9SBarry Smith   *m = a->m; *n = a->n;
111517ab2063SBarry Smith   return 0;
111617ab2063SBarry Smith }
111717ab2063SBarry Smith 
11185615d1e5SSatish Balay #undef __FUNC__
11195615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11208f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
112117ab2063SBarry Smith {
1122416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1123416022c9SBarry Smith   *m = 0; *n = a->m;
112417ab2063SBarry Smith   return 0;
112517ab2063SBarry Smith }
1126026e39d0SSatish Balay 
11275615d1e5SSatish Balay #undef __FUNC__
11285615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11294e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
113017ab2063SBarry Smith {
1131416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1132c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
113317ab2063SBarry Smith 
1134e3372554SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,0,"Row out of range");
113517ab2063SBarry Smith 
1136416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1137416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
113817ab2063SBarry Smith   if (idx) {
1139416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
11404e093b46SBarry Smith     if (*nz && shift) {
11410452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
114217ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
11434e093b46SBarry Smith     } else if (*nz) {
11444e093b46SBarry Smith       *idx = itmp;
114517ab2063SBarry Smith     }
114617ab2063SBarry Smith     else *idx = 0;
114717ab2063SBarry Smith   }
114817ab2063SBarry Smith   return 0;
114917ab2063SBarry Smith }
115017ab2063SBarry Smith 
11515615d1e5SSatish Balay #undef __FUNC__
11525615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
11534e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
115417ab2063SBarry Smith {
11554e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11564e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
115717ab2063SBarry Smith   return 0;
115817ab2063SBarry Smith }
115917ab2063SBarry Smith 
11605615d1e5SSatish Balay #undef __FUNC__
11615615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
11628f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
116317ab2063SBarry Smith {
1164416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1165416022c9SBarry Smith   Scalar     *v = a->a;
116617ab2063SBarry Smith   double     sum = 0.0;
1167416022c9SBarry Smith   int        i, j,shift = a->indexshift;
116817ab2063SBarry Smith 
116917ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1170416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
117117ab2063SBarry Smith #if defined(PETSC_COMPLEX)
117217ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
117317ab2063SBarry Smith #else
117417ab2063SBarry Smith       sum += (*v)*(*v); v++;
117517ab2063SBarry Smith #endif
117617ab2063SBarry Smith     }
117717ab2063SBarry Smith     *norm = sqrt(sum);
117817ab2063SBarry Smith   }
117917ab2063SBarry Smith   else if (type == NORM_1) {
118017ab2063SBarry Smith     double *tmp;
1181416022c9SBarry Smith     int    *jj = a->j;
118266963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1183cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
118417ab2063SBarry Smith     *norm = 0.0;
1185416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1186a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
118717ab2063SBarry Smith     }
1188416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
118917ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
119017ab2063SBarry Smith     }
11910452661fSBarry Smith     PetscFree(tmp);
119217ab2063SBarry Smith   }
119317ab2063SBarry Smith   else if (type == NORM_INFINITY) {
119417ab2063SBarry Smith     *norm = 0.0;
1195416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1196416022c9SBarry Smith       v = a->a + a->i[j] + shift;
119717ab2063SBarry Smith       sum = 0.0;
1198416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1199cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
120017ab2063SBarry Smith       }
120117ab2063SBarry Smith       if (sum > *norm) *norm = sum;
120217ab2063SBarry Smith     }
120317ab2063SBarry Smith   }
120417ab2063SBarry Smith   else {
1205e3372554SBarry Smith     SETERRQ(1,0,"No support for two norm yet");
120617ab2063SBarry Smith   }
120717ab2063SBarry Smith   return 0;
120817ab2063SBarry Smith }
120917ab2063SBarry Smith 
12105615d1e5SSatish Balay #undef __FUNC__
12115615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12128f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
121317ab2063SBarry Smith {
1214416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1215416022c9SBarry Smith   Mat        C;
1216416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1217416022c9SBarry Smith   int        shift = a->indexshift;
1218d5d45c9bSBarry Smith   Scalar     *array = a->a;
121917ab2063SBarry Smith 
12203638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
1221e3372554SBarry Smith     SETERRQ(1,0,"Square matrix only for in-place");
12220452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1223cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
122417ab2063SBarry Smith   if (shift) {
122517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
122617ab2063SBarry Smith   }
122717ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1228416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12290452661fSBarry Smith   PetscFree(col);
123017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
123117ab2063SBarry Smith     len = ai[i+1]-ai[i];
1232416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
123317ab2063SBarry Smith     array += len; aj += len;
123417ab2063SBarry Smith   }
123517ab2063SBarry Smith   if (shift) {
123617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
123717ab2063SBarry Smith   }
123817ab2063SBarry Smith 
12396d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12406d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
124117ab2063SBarry Smith 
12423638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1243416022c9SBarry Smith     *B = C;
124417ab2063SBarry Smith   } else {
1245416022c9SBarry Smith     /* This isn't really an in-place transpose */
12460452661fSBarry Smith     PetscFree(a->a);
12470452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
12480452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
12490452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
12500452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
12510452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
12521073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
12530452661fSBarry Smith     PetscFree(a);
1254f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
12550452661fSBarry Smith     PetscHeaderDestroy(C);
125617ab2063SBarry Smith   }
125717ab2063SBarry Smith   return 0;
125817ab2063SBarry Smith }
125917ab2063SBarry Smith 
12605615d1e5SSatish Balay #undef __FUNC__
12615615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
12628f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
126317ab2063SBarry Smith {
1264416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
126517ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1266d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
126717ab2063SBarry Smith 
126817ab2063SBarry Smith   if (ll) {
12693ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
12703ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
12719b1297e1SSatish Balay     VecGetLocalSize_Fast(ll,m);
1272e3372554SBarry Smith     if (m != a->m) SETERRQ(1,0,"Left scaling vector wrong length");
127390f02eecSBarry Smith     VecGetArray_Fast(ll,l);
1274416022c9SBarry Smith     v = a->a;
127517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
127617ab2063SBarry Smith       x = l[i];
1277416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
127817ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
127917ab2063SBarry Smith     }
128044cd7ae7SLois Curfman McInnes     PLogFlops(nz);
128117ab2063SBarry Smith   }
128217ab2063SBarry Smith   if (rr) {
12839b1297e1SSatish Balay     VecGetLocalSize_Fast(rr,n);
1284e3372554SBarry Smith     if (n != a->n) SETERRQ(1,0,"Right scaling vector wrong length");
128590f02eecSBarry Smith     VecGetArray_Fast(rr,r);
1286416022c9SBarry Smith     v = a->a; jj = a->j;
128717ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
128817ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
128917ab2063SBarry Smith     }
129044cd7ae7SLois Curfman McInnes     PLogFlops(nz);
129117ab2063SBarry Smith   }
129217ab2063SBarry Smith   return 0;
129317ab2063SBarry Smith }
129417ab2063SBarry Smith 
12955615d1e5SSatish Balay #undef __FUNC__
12965615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
12978f6be9afSLois Curfman McInnes int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
129817ab2063SBarry Smith {
1299db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
130002834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
130199141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1302a2744918SBarry Smith   register int sum,lensi;
130399141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
130499141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
130599141d43SSatish Balay   Scalar       *a_new,*mat_a;
1306416022c9SBarry Smith   Mat          C;
130717ab2063SBarry Smith 
1308b48a1e75SSatish Balay   ierr = ISSorted(isrow,(PetscTruth*)&i);
1309e3372554SBarry Smith   if (!i) SETERRQ(1,0,"ISrow is not sorted");
131099141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
1311e3372554SBarry Smith   if (!i) SETERRQ(1,0,"IScol is not sorted");
131299141d43SSatish Balay 
131317ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
131417ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
131517ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
131617ab2063SBarry Smith 
13177264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
131802834360SBarry Smith     /* special case of contiguous rows */
131957faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
132002834360SBarry Smith     starts = lens + ncols;
132102834360SBarry Smith     /* loop over new rows determining lens and starting points */
132202834360SBarry Smith     for (i=0; i<nrows; i++) {
1323a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1324a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
132502834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1326d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
132702834360SBarry Smith           starts[i] = k;
132802834360SBarry Smith           break;
132902834360SBarry Smith 	}
133002834360SBarry Smith       }
1331a2744918SBarry Smith       sum = 0;
133202834360SBarry Smith       while (k < kend) {
1333d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1334a2744918SBarry Smith         sum++;
133502834360SBarry Smith       }
1336a2744918SBarry Smith       lens[i] = sum;
133702834360SBarry Smith     }
133802834360SBarry Smith     /* create submatrix */
1339cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
134008480c60SBarry Smith       int n_cols,n_rows;
134108480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1342e3372554SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,0,"");
1343d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
134408480c60SBarry Smith       C = *B;
134508480c60SBarry Smith     }
134608480c60SBarry Smith     else {
134702834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
134808480c60SBarry Smith     }
1349db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1350db02288aSLois Curfman McInnes 
135102834360SBarry Smith     /* loop over rows inserting into submatrix */
1352db02288aSLois Curfman McInnes     a_new    = c->a;
1353db02288aSLois Curfman McInnes     j_new    = c->j;
1354db02288aSLois Curfman McInnes     i_new    = c->i;
1355db02288aSLois Curfman McInnes     i_new[0] = -shift;
135602834360SBarry Smith     for (i=0; i<nrows; i++) {
1357a2744918SBarry Smith       ii    = starts[i];
1358a2744918SBarry Smith       lensi = lens[i];
1359a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1360a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
136102834360SBarry Smith       }
1362a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1363a2744918SBarry Smith       a_new      += lensi;
1364a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1365a2744918SBarry Smith       c->ilen[i]  = lensi;
136602834360SBarry Smith     }
13670452661fSBarry Smith     PetscFree(lens);
136802834360SBarry Smith   }
136902834360SBarry Smith   else {
137002834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
13710452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
137202834360SBarry Smith     ssmap = smap + shift;
137399141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1374cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
137517ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
137602834360SBarry Smith     /* determine lens of each row */
137702834360SBarry Smith     for (i=0; i<nrows; i++) {
1378d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
137902834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
138002834360SBarry Smith       lens[i] = 0;
138102834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1382d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
138302834360SBarry Smith           lens[i]++;
138402834360SBarry Smith         }
138502834360SBarry Smith       }
138602834360SBarry Smith     }
138717ab2063SBarry Smith     /* Create and fill new matrix */
1388a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
138999141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
139099141d43SSatish Balay 
1391e3372554SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(1,0,"");
139299141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1393e3372554SBarry Smith         SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros");
139499141d43SSatish Balay       }
139599141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
139608480c60SBarry Smith       C = *B;
139708480c60SBarry Smith     }
139808480c60SBarry Smith     else {
139902834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
140008480c60SBarry Smith     }
140199141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
140217ab2063SBarry Smith     for (i=0; i<nrows; i++) {
140399141d43SSatish Balay       row    = irow[i];
140417ab2063SBarry Smith       nznew  = 0;
140599141d43SSatish Balay       kstart = ai[row]+shift;
140699141d43SSatish Balay       kend   = kstart + a->ilen[row];
140799141d43SSatish Balay       mat_i  = c->i[i]+shift;
140899141d43SSatish Balay       mat_j  = c->j + mat_i;
140999141d43SSatish Balay       mat_a  = c->a + mat_i;
141099141d43SSatish Balay       mat_ilen = c->ilen + i;
141117ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
141299141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
141399141d43SSatish Balay           *mat_j++ = tcol - (!shift);
141499141d43SSatish Balay           *mat_a++ = a->a[k];
141599141d43SSatish Balay           (*mat_ilen)++;
141699141d43SSatish Balay 
141717ab2063SBarry Smith         }
141817ab2063SBarry Smith       }
141917ab2063SBarry Smith     }
142002834360SBarry Smith     /* Free work space */
142102834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
142299141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
142302834360SBarry Smith   }
14246d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
14256d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
142617ab2063SBarry Smith 
142717ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1428416022c9SBarry Smith   *B = C;
142917ab2063SBarry Smith   return 0;
143017ab2063SBarry Smith }
143117ab2063SBarry Smith 
1432a871dcd8SBarry Smith /*
143363b91edcSBarry Smith      note: This can only work for identity for row and col. It would
143463b91edcSBarry Smith    be good to check this and otherwise generate an error.
1435a871dcd8SBarry Smith */
14365615d1e5SSatish Balay #undef __FUNC__
14375615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
14388f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1439a871dcd8SBarry Smith {
144063b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
144108480c60SBarry Smith   int        ierr;
144263b91edcSBarry Smith   Mat        outA;
144363b91edcSBarry Smith 
1444e3372554SBarry Smith   if (fill != 0) SETERRQ(1,0,"Only fill=0 supported");
1445a871dcd8SBarry Smith 
144663b91edcSBarry Smith   outA          = inA;
144763b91edcSBarry Smith   inA->factor   = FACTOR_LU;
144863b91edcSBarry Smith   a->row        = row;
144963b91edcSBarry Smith   a->col        = col;
145063b91edcSBarry Smith 
145194a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
14520452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
145394a9d846SBarry Smith   }
145463b91edcSBarry Smith 
145508480c60SBarry Smith   if (!a->diag) {
145608480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
145763b91edcSBarry Smith   }
145863b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1459a871dcd8SBarry Smith   return 0;
1460a871dcd8SBarry Smith }
1461a871dcd8SBarry Smith 
1462f0b747eeSBarry Smith #include "pinclude/plapack.h"
14635615d1e5SSatish Balay #undef __FUNC__
14645615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
14658f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1466f0b747eeSBarry Smith {
1467f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1468f0b747eeSBarry Smith   int        one = 1;
1469f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1470f0b747eeSBarry Smith   PLogFlops(a->nz);
1471f0b747eeSBarry Smith   return 0;
1472f0b747eeSBarry Smith }
1473f0b747eeSBarry Smith 
14745615d1e5SSatish Balay #undef __FUNC__
14755615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
14768f6be9afSLois Curfman McInnes int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1477cddf8d76SBarry Smith                                     Mat **B)
1478cddf8d76SBarry Smith {
1479cddf8d76SBarry Smith   int ierr,i;
1480cddf8d76SBarry Smith 
1481cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
14820452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1483cddf8d76SBarry Smith   }
1484cddf8d76SBarry Smith 
1485cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1486905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1487cddf8d76SBarry Smith   }
1488cddf8d76SBarry Smith   return 0;
1489cddf8d76SBarry Smith }
1490cddf8d76SBarry Smith 
14915615d1e5SSatish Balay #undef __FUNC__
14925615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
14938f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
14945a838052SSatish Balay {
14955a838052SSatish Balay   *bs = 1;
14965a838052SSatish Balay   return 0;
14975a838052SSatish Balay }
14985a838052SSatish Balay 
14995615d1e5SSatish Balay #undef __FUNC__
15005615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15018f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15024dcbc457SBarry Smith {
1503e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
150406763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15058a047759SSatish Balay   int        start, end, *ai, *aj;
150606763907SSatish Balay   char       *table;
15078a047759SSatish Balay   shift = a->indexshift;
1508e4d965acSSatish Balay   m     = a->m;
1509e4d965acSSatish Balay   ai    = a->i;
15108a047759SSatish Balay   aj    = a->j+shift;
15118a047759SSatish Balay 
1512e3372554SBarry Smith   if (ov < 0)  SETERRQ(1,0,"illegal overlap value used");
151306763907SSatish Balay 
151406763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
151506763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
151606763907SSatish Balay 
1517e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1518b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1519e4d965acSSatish Balay     isz  = 0;
152006763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1521e4d965acSSatish Balay 
1522e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
15234dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
152477c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1525e4d965acSSatish Balay 
1526dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1527e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
152806763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
15294dcbc457SBarry Smith     }
153006763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
153106763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1532e4d965acSSatish Balay 
153304a348a9SBarry Smith     k = 0;
153404a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
153504a348a9SBarry Smith       n = isz;
153606763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1537e4d965acSSatish Balay         row   = nidx[k];
1538e4d965acSSatish Balay         start = ai[row];
1539e4d965acSSatish Balay         end   = ai[row+1];
154004a348a9SBarry Smith         for ( l = start; l<end ; l++){
15418a047759SSatish Balay           val = aj[l] + shift;
154206763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1543e4d965acSSatish Balay         }
1544e4d965acSSatish Balay       }
1545e4d965acSSatish Balay     }
1546029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1547e4d965acSSatish Balay   }
154804a348a9SBarry Smith   PetscFree(table);
154906763907SSatish Balay   PetscFree(nidx);
1550e4d965acSSatish Balay   return 0;
15514dcbc457SBarry Smith }
155217ab2063SBarry Smith 
15530513a670SBarry Smith /* -------------------------------------------------------------- */
15540513a670SBarry Smith #undef __FUNC__
15550513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
15560513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
15570513a670SBarry Smith {
15580513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
15590513a670SBarry Smith   Scalar     *vwork;
15600513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
15610513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
15620513a670SBarry Smith 
15630513a670SBarry Smith   ierr = ISGetIndices(rowp,&row); CHKERRQ(ierr);
15640513a670SBarry Smith   ierr = ISGetIndices(colp,&col); CHKERRQ(ierr);
15650513a670SBarry Smith 
15660513a670SBarry Smith   /* determine lengths of permuted rows */
15670513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
15680513a670SBarry Smith   for (i=0; i<m; i++ ) {
15690513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
15700513a670SBarry Smith   }
15710513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
15720513a670SBarry Smith   PetscFree(lens);
15730513a670SBarry Smith 
15740513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
15750513a670SBarry Smith   for (i=0; i<m; i++) {
15760513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15770513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
15780513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
15790513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15800513a670SBarry Smith   }
15810513a670SBarry Smith   PetscFree(cnew);
15820513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
15830513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
15840513a670SBarry Smith   ierr = ISRestoreIndices(rowp,&row); CHKERRQ(ierr);
15850513a670SBarry Smith   ierr = ISRestoreIndices(colp,&col); CHKERRQ(ierr);
15860513a670SBarry Smith   return 0;
15870513a670SBarry Smith }
15880513a670SBarry Smith 
15895615d1e5SSatish Balay #undef __FUNC__
15905615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1591682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1592682d7d0cSBarry Smith {
1593682d7d0cSBarry Smith   static int called = 0;
1594682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1595682d7d0cSBarry Smith 
1596682d7d0cSBarry Smith   if (called) return 0; else called = 1;
159777c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
15980f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
15990f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
16000f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_no_inode: Do not use inodes\n");
16010f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1602682d7d0cSBarry Smith #if defined(HAVE_ESSL)
16030f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1604682d7d0cSBarry Smith #endif
1605682d7d0cSBarry Smith   return 0;
1606682d7d0cSBarry Smith }
16078f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1608a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1609a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1610a93ec695SBarry Smith 
1611682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
161217ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
161317ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1614416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1615416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
161617ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
161717ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
161817ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
161917ab2063SBarry Smith        MatRelax_SeqAIJ,
162017ab2063SBarry Smith        MatTranspose_SeqAIJ,
16217264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1622f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
162317ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
162417ab2063SBarry Smith        MatCompress_SeqAIJ,
162517ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
162617ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
162717ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
162817ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
162994a9d846SBarry Smith        0,0,
16303d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1631cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
16327eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1633682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1634f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
16355a838052SSatish Balay        MatScale_SeqAIJ,0,0,
16366945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
16376945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
16383b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
16393b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
16403b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1641a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1642a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
16430513a670SBarry Smith        MatColoringPatch_SeqAIJ,
16440513a670SBarry Smith        0,
16450513a670SBarry Smith        MatPermute_SeqAIJ};
164617ab2063SBarry Smith 
164717ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
164817ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
164917ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
165017ab2063SBarry Smith 
16515615d1e5SSatish Balay #undef __FUNC__
16525615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
165317ab2063SBarry Smith /*@C
1654682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
16550d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
16566e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
16572bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
16582bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
165917ab2063SBarry Smith 
166017ab2063SBarry Smith    Input Parameters:
1661029af93fSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
166217ab2063SBarry Smith .  m - number of rows
166317ab2063SBarry Smith .  n - number of columns
166417ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
16652bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
16662bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
166717ab2063SBarry Smith 
166817ab2063SBarry Smith    Output Parameter:
1669416022c9SBarry Smith .  A - the matrix
167017ab2063SBarry Smith 
167117ab2063SBarry Smith    Notes:
167217ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
167317ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
16740002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
167544cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
167617ab2063SBarry Smith 
167717ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1678a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
16793d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
16806da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
168117ab2063SBarry Smith 
1682682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1683682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1684682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
16856c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
16866c7ebb05SLois Curfman McInnes 
16876c7ebb05SLois Curfman McInnes    Options Database Keys:
16886c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
16896e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
16906e62573dSLois Curfman McInnes $        (max limit=5)
16916e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
16926e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
16936e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
169417ab2063SBarry Smith 
169517ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
169617ab2063SBarry Smith @*/
1697416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
169817ab2063SBarry Smith {
1699416022c9SBarry Smith   Mat        B;
1700416022c9SBarry Smith   Mat_SeqAIJ *b;
17016945ee14SBarry Smith   int        i, len, ierr, flg,size;
17026945ee14SBarry Smith 
17036945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1704e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
1705d5d45c9bSBarry Smith 
1706416022c9SBarry Smith   *A                  = 0;
1707f09e8eb9SSatish Balay   PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1708416022c9SBarry Smith   PLogObjectCreate(B);
17090452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
171044cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1711416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1712416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1713416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1714416022c9SBarry Smith   B->factor           = 0;
1715416022c9SBarry Smith   B->lupivotthreshold = 1.0;
171690f02eecSBarry Smith   B->mapping          = 0;
17177a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
171869957df2SSatish Balay                           &flg); CHKERRQ(ierr);
17197a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
17207a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
17217a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1722416022c9SBarry Smith   b->row              = 0;
1723416022c9SBarry Smith   b->col              = 0;
1724416022c9SBarry Smith   b->indexshift       = 0;
1725b810aeb4SBarry Smith   b->reallocs         = 0;
172669957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
172769957df2SSatish Balay   if (flg) b->indexshift = -1;
172817ab2063SBarry Smith 
172944cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
173044cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
17310452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1732b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
17337b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
17347b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1735416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
173617ab2063SBarry Smith     nz = nz*m;
173717ab2063SBarry Smith   }
173817ab2063SBarry Smith   else {
173917ab2063SBarry Smith     nz = 0;
1740416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
174117ab2063SBarry Smith   }
174217ab2063SBarry Smith 
174317ab2063SBarry Smith   /* allocate the matrix space */
1744416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
17450452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1746416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1747cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1748416022c9SBarry Smith   b->i  = b->j + nz;
1749416022c9SBarry Smith   b->singlemalloc = 1;
175017ab2063SBarry Smith 
1751416022c9SBarry Smith   b->i[0] = -b->indexshift;
175217ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1753416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
175417ab2063SBarry Smith   }
175517ab2063SBarry Smith 
1756416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
17570452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1758f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1759416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
176017ab2063SBarry Smith 
1761416022c9SBarry Smith   b->nz               = 0;
1762416022c9SBarry Smith   b->maxnz            = nz;
1763416022c9SBarry Smith   b->sorted           = 0;
1764416022c9SBarry Smith   b->roworiented      = 1;
1765416022c9SBarry Smith   b->nonew            = 0;
1766416022c9SBarry Smith   b->diag             = 0;
1767416022c9SBarry Smith   b->solve_work       = 0;
176871bd300dSLois Curfman McInnes   b->spptr            = 0;
1769754ec7b1SSatish Balay   b->inode.node_count = 0;
1770754ec7b1SSatish Balay   b->inode.size       = 0;
17716c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
17726c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
17734e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
177417ab2063SBarry Smith 
1775416022c9SBarry Smith   *A = B;
17764e220ebcSLois Curfman McInnes 
17774b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
17784b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
177969957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
178069957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
17814b14c69eSBarry Smith #endif
178269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
178369957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
178469957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
178569957df2SSatish Balay   if (flg) {
1786e3372554SBarry Smith     if (!b->indexshift) SETERRQ(1,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1787416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
178817ab2063SBarry Smith   }
178969957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
179069957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
179117ab2063SBarry Smith   return 0;
179217ab2063SBarry Smith }
179317ab2063SBarry Smith 
17945615d1e5SSatish Balay #undef __FUNC__
17955615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
17963d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
179717ab2063SBarry Smith {
1798416022c9SBarry Smith   Mat        C;
1799416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
180008480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
180117ab2063SBarry Smith 
18024043dd9cSLois Curfman McInnes   *B = 0;
1803f09e8eb9SSatish Balay   PetscHeaderCreate(C,_p_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1804416022c9SBarry Smith   PLogObjectCreate(C);
18050452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
180641c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1807416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1808416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1809416022c9SBarry Smith   C->factor     = A->factor;
1810416022c9SBarry Smith   c->row        = 0;
1811416022c9SBarry Smith   c->col        = 0;
1812416022c9SBarry Smith   c->indexshift = shift;
1813c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
181417ab2063SBarry Smith 
181544cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
181644cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
181744cd7ae7SLois Curfman McInnes   C->M          = a->m;
181844cd7ae7SLois Curfman McInnes   C->N          = a->n;
181917ab2063SBarry Smith 
18200452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
18210452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
182217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1823416022c9SBarry Smith     c->imax[i] = a->imax[i];
1824416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
182517ab2063SBarry Smith   }
182617ab2063SBarry Smith 
182717ab2063SBarry Smith   /* allocate the matrix space */
1828416022c9SBarry Smith   c->singlemalloc = 1;
1829416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
18300452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1831416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1832416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1833416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
183417ab2063SBarry Smith   if (m > 0) {
1835416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
183608480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1837416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
183817ab2063SBarry Smith     }
183908480c60SBarry Smith   }
184017ab2063SBarry Smith 
1841f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1842416022c9SBarry Smith   c->sorted      = a->sorted;
1843416022c9SBarry Smith   c->roworiented = a->roworiented;
1844416022c9SBarry Smith   c->nonew       = a->nonew;
18457a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
184617ab2063SBarry Smith 
1847416022c9SBarry Smith   if (a->diag) {
18480452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1849416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
185017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1851416022c9SBarry Smith       c->diag[i] = a->diag[i];
185217ab2063SBarry Smith     }
185317ab2063SBarry Smith   }
1854416022c9SBarry Smith   else c->diag          = 0;
18556c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
18566c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1857754ec7b1SSatish Balay   if (a->inode.size){
1858daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1859754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1860daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1861754ec7b1SSatish Balay   } else {
1862754ec7b1SSatish Balay     c->inode.size       = 0;
1863754ec7b1SSatish Balay     c->inode.node_count = 0;
1864754ec7b1SSatish Balay   }
1865416022c9SBarry Smith   c->nz                 = a->nz;
1866416022c9SBarry Smith   c->maxnz              = a->maxnz;
1867416022c9SBarry Smith   c->solve_work         = 0;
186876dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1869754ec7b1SSatish Balay 
1870416022c9SBarry Smith   *B = C;
187117ab2063SBarry Smith   return 0;
187217ab2063SBarry Smith }
187317ab2063SBarry Smith 
18745615d1e5SSatish Balay #undef __FUNC__
18755615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
187619bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
187717ab2063SBarry Smith {
1878416022c9SBarry Smith   Mat_SeqAIJ   *a;
1879416022c9SBarry Smith   Mat          B;
188017699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1881bcd2baecSBarry Smith   MPI_Comm     comm;
188217ab2063SBarry Smith 
188319bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
188417699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
1885e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"view must have one processor");
188690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
188777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1888e3372554SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object in file");
188917ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
189017ab2063SBarry Smith 
189117ab2063SBarry Smith   /* read in row lengths */
18920452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
189377c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
189417ab2063SBarry Smith 
189517ab2063SBarry Smith   /* create our matrix */
1896416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1897416022c9SBarry Smith   B = *A;
1898416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1899416022c9SBarry Smith   shift = a->indexshift;
190017ab2063SBarry Smith 
190117ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
190277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
190317ab2063SBarry Smith   if (shift) {
190417ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1905416022c9SBarry Smith       a->j[i] += 1;
190617ab2063SBarry Smith     }
190717ab2063SBarry Smith   }
190817ab2063SBarry Smith 
190917ab2063SBarry Smith   /* read in nonzero values */
191077c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
191117ab2063SBarry Smith 
191217ab2063SBarry Smith   /* set matrix "i" values */
1913416022c9SBarry Smith   a->i[0] = -shift;
191417ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1915416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1916416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
191717ab2063SBarry Smith   }
19180452661fSBarry Smith   PetscFree(rowlengths);
191917ab2063SBarry Smith 
19206d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
19216d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
192217ab2063SBarry Smith   return 0;
192317ab2063SBarry Smith }
192417ab2063SBarry Smith 
19255615d1e5SSatish Balay #undef __FUNC__
19265615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
19278f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
19287264ac53SSatish Balay {
19297264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
19307264ac53SSatish Balay 
1931e3372554SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,0,"Matrices must be same type");
19327264ac53SSatish Balay 
19337264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
19347264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1935bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
193677c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1937bcd2baecSBarry Smith   }
19387264ac53SSatish Balay 
19397264ac53SSatish Balay   /* if the a->i are the same */
1940bcd2baecSBarry Smith   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
194177c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
19427264ac53SSatish Balay   }
19437264ac53SSatish Balay 
19447264ac53SSatish Balay   /* if a->j are the same */
1945bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
194677c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1947bcd2baecSBarry Smith   }
1948bcd2baecSBarry Smith 
1949bcd2baecSBarry Smith   /* if a->a are the same */
195019bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
195177c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
19527264ac53SSatish Balay   }
195377c4ece6SBarry Smith   *flg = PETSC_TRUE;
19547264ac53SSatish Balay   return 0;
19557264ac53SSatish Balay 
19567264ac53SSatish Balay }
1957