xref: /petsc/src/mat/impls/aij/seq/aij.c (revision c2653b3d308917c4b5698e5f59fd76add26d7f41)
16d84be18SBarry Smith 
26945ee14SBarry Smith 
317ab2063SBarry Smith #ifndef lint
4*c2653b3dSLois Curfman McInnes static char vcid[] = "$Id: aij.c,v 1.209 1997/03/01 15:48:57 bsmith Exp bsmith $";
517ab2063SBarry Smith #endif
617ab2063SBarry Smith 
7d5d45c9bSBarry Smith /*
85a838052SSatish Balay B    Defines the basic matrix operations for the AIJ (compressed row)
9d5d45c9bSBarry Smith   matrix storage format.
10d5d45c9bSBarry Smith */
1170f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
12f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
13f5eb4b81SSatish Balay #include "src/inline/spops.h"
14e4d965acSSatish Balay #include "petsc.h"
15f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1617ab2063SBarry Smith 
17a2ce50c7SBarry Smith /*
18a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
19a2ce50c7SBarry Smith */
205615d1e5SSatish Balay #undef __FUNC__
215615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ"
22a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
23a2ce50c7SBarry Smith {
24a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
25a2ce50c7SBarry Smith   int        ierr = 1;
26a2ce50c7SBarry Smith 
27e3372554SBarry Smith   SETERRQ(ierr,0,"Not implemented");
28a2ce50c7SBarry Smith }
29a2ce50c7SBarry Smith 
30bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3117ab2063SBarry Smith 
325615d1e5SSatish Balay #undef __FUNC__
335615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ"
3431625ec5SSatish Balay static int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
356945ee14SBarry Smith                            PetscTruth *done)
3617ab2063SBarry Smith {
37416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
386945ee14SBarry Smith   int        ierr,i,ishift;
3917ab2063SBarry Smith 
4031625ec5SSatish Balay   *m     = A->m;
416945ee14SBarry Smith   if (!ia) return 0;
426945ee14SBarry Smith   ishift = a->indexshift;
436945ee14SBarry Smith   if (symmetric) {
4431625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
456945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
4631625ec5SSatish Balay     int nz = a->i[a->m];
473b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
4831625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
493b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
503b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
5131625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
526945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
5331625ec5SSatish Balay     int nz = a->i[a->m] + 1;
543b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
5531625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
563b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
573b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
5831625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
596945ee14SBarry Smith   } else {
606945ee14SBarry Smith     *ia = a->i; *ja = a->j;
61a2ce50c7SBarry Smith   }
62a2ce50c7SBarry Smith 
63a2744918SBarry Smith   return 0;
64a2744918SBarry Smith }
65a2744918SBarry Smith 
665615d1e5SSatish Balay #undef __FUNC__
675615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ"
683b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
696945ee14SBarry Smith                                PetscTruth *done)
706945ee14SBarry Smith {
716945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
723b2fbd54SBarry Smith   int        ishift = a->indexshift;
736945ee14SBarry Smith 
746945ee14SBarry Smith   if (!ia) return 0;
753b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
766945ee14SBarry Smith     PetscFree(*ia);
776945ee14SBarry Smith     PetscFree(*ja);
78bcd2baecSBarry Smith   }
7917ab2063SBarry Smith   return 0;
8017ab2063SBarry Smith }
8117ab2063SBarry Smith 
825615d1e5SSatish Balay #undef __FUNC__
835615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ"
8443a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
853b2fbd54SBarry Smith                            PetscTruth *done)
863b2fbd54SBarry Smith {
873b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
88a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
89a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
903b2fbd54SBarry Smith 
913b2fbd54SBarry Smith   *nn     = A->n;
923b2fbd54SBarry Smith   if (!ia) return 0;
933b2fbd54SBarry Smith   if (symmetric) {
94179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
953b2fbd54SBarry Smith   } else {
9661d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
973b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
983b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
99a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
1003b2fbd54SBarry Smith     jj = a->j;
1013b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
1023b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
1033b2fbd54SBarry Smith     }
1043b2fbd54SBarry Smith     cia[0] = oshift;
1053b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
1063b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
1073b2fbd54SBarry Smith     }
1083b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1093b2fbd54SBarry Smith     jj = a->j;
110a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
111a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
112a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1133b2fbd54SBarry Smith         col = *jj++ + ishift;
1143b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1153b2fbd54SBarry Smith       }
1163b2fbd54SBarry Smith     }
1173b2fbd54SBarry Smith     PetscFree(collengths);
1183b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1193b2fbd54SBarry Smith   }
1203b2fbd54SBarry Smith 
1213b2fbd54SBarry Smith   return 0;
1223b2fbd54SBarry Smith }
1233b2fbd54SBarry Smith 
1245615d1e5SSatish Balay #undef __FUNC__
1255615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ"
12643a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1273b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1283b2fbd54SBarry Smith {
1293b2fbd54SBarry Smith   if (!ia) return 0;
1303b2fbd54SBarry Smith 
1313b2fbd54SBarry Smith   PetscFree(*ia);
1323b2fbd54SBarry Smith   PetscFree(*ja);
1333b2fbd54SBarry Smith 
1343b2fbd54SBarry Smith   return 0;
1353b2fbd54SBarry Smith }
1363b2fbd54SBarry Smith 
137227d817aSBarry Smith #define CHUNKSIZE   15
13817ab2063SBarry Smith 
1395615d1e5SSatish Balay #undef __FUNC__
1405615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ"
14161d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
14217ab2063SBarry Smith {
143416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
144416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1454b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
146d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
147416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
14817ab2063SBarry Smith 
14917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
150416022c9SBarry Smith     row  = im[k];
1513b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
152e3372554SBarry Smith     if (row < 0) SETERRQ(1,0,"Negative row");
153e3372554SBarry Smith     if (row >= a->m) SETERRQ(1,0,"Row too large");
1543b2fbd54SBarry Smith #endif
15517ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
15617ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
157416022c9SBarry Smith     low = 0;
15817ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1593b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
160e3372554SBarry Smith       if (in[l] < 0) SETERRQ(1,0,"Negative column");
161e3372554SBarry Smith       if (in[l] >= a->n) SETERRQ(1,0,"Column too large");
1623b2fbd54SBarry Smith #endif
1634b0e389bSBarry Smith       col = in[l] - shift;
1644b0e389bSBarry Smith       if (roworiented) {
1654b0e389bSBarry Smith         value = *v++;
1664b0e389bSBarry Smith       }
1674b0e389bSBarry Smith       else {
1684b0e389bSBarry Smith         value = v[k + l*m];
1694b0e389bSBarry Smith       }
170416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
171416022c9SBarry Smith       while (high-low > 5) {
172416022c9SBarry Smith         t = (low+high)/2;
173416022c9SBarry Smith         if (rp[t] > col) high = t;
174416022c9SBarry Smith         else             low  = t;
17517ab2063SBarry Smith       }
176416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
17717ab2063SBarry Smith         if (rp[i] > col) break;
17817ab2063SBarry Smith         if (rp[i] == col) {
179416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
18017ab2063SBarry Smith           else                  ap[i] = value;
18117ab2063SBarry Smith           goto noinsert;
18217ab2063SBarry Smith         }
18317ab2063SBarry Smith       }
184*c2653b3dSLois Curfman McInnes       if (nonew == 1) goto noinsert;
185*c2653b3dSLois Curfman McInnes       else if (nonew == -1) SETERRQ(1,1,"Inserting a new nonzero in the matrix");
18617ab2063SBarry Smith       if (nrow >= rmax) {
18717ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
188416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
18917ab2063SBarry Smith         Scalar *new_a;
19017ab2063SBarry Smith 
19117ab2063SBarry Smith         /* malloc new storage space */
192416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1930452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
19417ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
19517ab2063SBarry Smith         new_i   = new_j + new_nz;
19617ab2063SBarry Smith 
19717ab2063SBarry Smith         /* copy over old data into new slots */
19817ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
199416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
200416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
201416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
202416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
20317ab2063SBarry Smith                                                            len*sizeof(int));
204416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
205416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
20617ab2063SBarry Smith                                                            len*sizeof(Scalar));
20717ab2063SBarry Smith         /* free up old matrix storage */
2080452661fSBarry Smith         PetscFree(a->a);
2090452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
210416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
211416022c9SBarry Smith         a->singlemalloc = 1;
21217ab2063SBarry Smith 
21317ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
214416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
215416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
216416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
217b810aeb4SBarry Smith         a->reallocs++;
21817ab2063SBarry Smith       }
219416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
220416022c9SBarry Smith       /* shift up all the later entries in this row */
221416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
22217ab2063SBarry Smith         rp[ii+1] = rp[ii];
22317ab2063SBarry Smith         ap[ii+1] = ap[ii];
22417ab2063SBarry Smith       }
22517ab2063SBarry Smith       rp[i] = col;
22617ab2063SBarry Smith       ap[i] = value;
22717ab2063SBarry Smith       noinsert:;
228416022c9SBarry Smith       low = i + 1;
22917ab2063SBarry Smith     }
23017ab2063SBarry Smith     ailen[row] = nrow;
23117ab2063SBarry Smith   }
23217ab2063SBarry Smith   return 0;
23317ab2063SBarry Smith }
23417ab2063SBarry Smith 
2355615d1e5SSatish Balay #undef __FUNC__
2365615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ"
2377eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2387eb43aa7SLois Curfman McInnes {
2397eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
240b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2417eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2427eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2437eb43aa7SLois Curfman McInnes 
2447eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2457eb43aa7SLois Curfman McInnes     row  = im[k];
246e3372554SBarry Smith     if (row < 0) SETERRQ(1,0,"Negative row");
247e3372554SBarry Smith     if (row >= a->m) SETERRQ(1,0,"Row too large");
2487eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2497eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2507eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
251e3372554SBarry Smith       if (in[l] < 0) SETERRQ(1,0,"Negative column");
252e3372554SBarry Smith       if (in[l] >= a->n) SETERRQ(1,0,"Column too large");
2537eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2547eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2557eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2567eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2577eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2587eb43aa7SLois Curfman McInnes         else             low  = t;
2597eb43aa7SLois Curfman McInnes       }
2607eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2617eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2627eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
263b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2647eb43aa7SLois Curfman McInnes           goto finished;
2657eb43aa7SLois Curfman McInnes         }
2667eb43aa7SLois Curfman McInnes       }
267b49de8d1SLois Curfman McInnes       *v++ = zero;
2687eb43aa7SLois Curfman McInnes       finished:;
2697eb43aa7SLois Curfman McInnes     }
2707eb43aa7SLois Curfman McInnes   }
2717eb43aa7SLois Curfman McInnes   return 0;
2727eb43aa7SLois Curfman McInnes }
2737eb43aa7SLois Curfman McInnes 
27417ab2063SBarry Smith #include "draw.h"
27517ab2063SBarry Smith #include "pinclude/pviewer.h"
27677c4ece6SBarry Smith #include "sys.h"
27717ab2063SBarry Smith 
2785615d1e5SSatish Balay #undef __FUNC__
2795615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary"
280416022c9SBarry Smith static 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"
315416022c9SBarry Smith static 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)
35844cd7ae7SLois 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]));
36044cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
36144cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
36244cd7ae7SLois Curfman McInnes #else
36344cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
36444cd7ae7SLois Curfman McInnes #endif
36544cd7ae7SLois Curfman McInnes       }
36644cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
36744cd7ae7SLois Curfman McInnes     }
36844cd7ae7SLois Curfman McInnes   }
36917ab2063SBarry Smith   else {
37017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
37117ab2063SBarry Smith       fprintf(fd,"row %d:",i);
372416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
37317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
374416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
375416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
37617ab2063SBarry Smith         }
37717ab2063SBarry Smith         else {
378416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
37917ab2063SBarry Smith         }
38017ab2063SBarry Smith #else
381416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
38217ab2063SBarry Smith #endif
38317ab2063SBarry Smith       }
38417ab2063SBarry Smith       fprintf(fd,"\n");
38517ab2063SBarry Smith     }
38617ab2063SBarry Smith   }
38717ab2063SBarry Smith   fflush(fd);
388416022c9SBarry Smith   return 0;
389416022c9SBarry Smith }
390416022c9SBarry Smith 
3915615d1e5SSatish Balay #undef __FUNC__
3925615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Draw"
393416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
394416022c9SBarry Smith {
395416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
396cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
3970513a670SBarry Smith   int         format;
39894a9d846SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r,maxv = 0.0;
399bcd2baecSBarry Smith   Draw        draw;
400cddf8d76SBarry Smith   DrawButton  button;
40119bcc07fSBarry Smith   PetscTruth  isnull;
402cddf8d76SBarry Smith 
4030513a670SBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
4040513a670SBarry Smith   ierr = DrawClear(draw); CHKERRQ(ierr);
4050513a670SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
40619bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
40719bcc07fSBarry Smith 
408416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
409416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
410416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
411416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
4120513a670SBarry Smith 
4130513a670SBarry Smith   if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4140513a670SBarry Smith     /* Blue for negative, Cyan for zero and  Red for positive */
415cddf8d76SBarry Smith     color = DRAW_BLUE;
416416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
417cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
418416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
419cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
420cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
421cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
422cddf8d76SBarry Smith #else
423cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
424cddf8d76SBarry Smith #endif
425cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
426cddf8d76SBarry Smith       }
427cddf8d76SBarry Smith     }
428cddf8d76SBarry Smith     color = DRAW_CYAN;
429cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
430cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
431cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
432cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
433cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
434cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
435cddf8d76SBarry Smith       }
436cddf8d76SBarry Smith     }
437cddf8d76SBarry Smith     color = DRAW_RED;
438cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
439cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
440cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
441cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
442cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
443cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
444cddf8d76SBarry Smith #else
445cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
446cddf8d76SBarry Smith #endif
447cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
448416022c9SBarry Smith       }
449416022c9SBarry Smith     }
4500513a670SBarry Smith   } else {
4510513a670SBarry Smith     /* use contour shading to indicate magnitude of values */
4520513a670SBarry Smith     /* first determine max of all nonzero values */
4530513a670SBarry Smith     int    nz = a->nz,count;
4540513a670SBarry Smith     Draw   popup;
4550513a670SBarry Smith 
4560513a670SBarry Smith     for ( i=0; i<nz; i++ ) {
4570513a670SBarry Smith       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
4580513a670SBarry Smith     }
4590513a670SBarry Smith     ierr = DrawCreatePopUp(draw,&popup); CHKERRQ(ierr);
4600513a670SBarry Smith     ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr);
4610513a670SBarry Smith     count = 0;
4620513a670SBarry Smith     for ( i=0; i<m; i++ ) {
4630513a670SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
4640513a670SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
4650513a670SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
4660513a670SBarry Smith         color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
4670513a670SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
4680513a670SBarry Smith         count++;
4690513a670SBarry Smith       }
4700513a670SBarry Smith     }
4710513a670SBarry Smith   }
472416022c9SBarry Smith   DrawFlush(draw);
473cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
474cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
475cddf8d76SBarry Smith 
476cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
4776945ee14SBarry Smith   ierr = DrawCheckResizedWindow(draw);
478cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
479cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
480cddf8d76SBarry Smith     DrawClear(draw);
481cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
482cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
483cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
484cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
485cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
486cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
487cddf8d76SBarry Smith     w *= scale; h *= scale;
488cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
4890513a670SBarry Smith     if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4900513a670SBarry Smith       /* Blue for negative, Cyan for zero and  Red for positive */
491cddf8d76SBarry Smith       color = DRAW_BLUE;
492cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
493cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
494cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
495cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
496cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
497cddf8d76SBarry Smith           if (real(a->a[j]) >=  0.) continue;
498cddf8d76SBarry Smith #else
499cddf8d76SBarry Smith           if (a->a[j] >=  0.) continue;
500cddf8d76SBarry Smith #endif
501cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
502cddf8d76SBarry Smith         }
503cddf8d76SBarry Smith       }
504cddf8d76SBarry Smith       color = DRAW_CYAN;
505cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
506cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
507cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
508cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
509cddf8d76SBarry Smith           if (a->a[j] !=  0.) continue;
510cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
511cddf8d76SBarry Smith         }
512cddf8d76SBarry Smith       }
513cddf8d76SBarry Smith       color = DRAW_RED;
514cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
515cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
516cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
517cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
518cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
519cddf8d76SBarry Smith           if (real(a->a[j]) <=  0.) continue;
520cddf8d76SBarry Smith #else
521cddf8d76SBarry Smith           if (a->a[j] <=  0.) continue;
522cddf8d76SBarry Smith #endif
523cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
524cddf8d76SBarry Smith         }
525cddf8d76SBarry Smith       }
5260513a670SBarry Smith     } else {
5270513a670SBarry Smith       /* use contour shading to indicate magnitude of values */
5280513a670SBarry Smith       int count = 0;
5290513a670SBarry Smith       for ( i=0; i<m; i++ ) {
5300513a670SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
5310513a670SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5320513a670SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
5330513a670SBarry Smith           color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5340513a670SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5350513a670SBarry Smith           count++;
5360513a670SBarry Smith         }
5370513a670SBarry Smith       }
5380513a670SBarry Smith     }
5390513a670SBarry Smith 
5406945ee14SBarry Smith     ierr = DrawCheckResizedWindow(draw);
541cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
542cddf8d76SBarry Smith   }
543416022c9SBarry Smith   return 0;
544416022c9SBarry Smith }
545416022c9SBarry Smith 
5465615d1e5SSatish Balay #undef __FUNC__
547c22c1629SBarry Smith #define __FUNC__ "MatView_SeqAIJ" /* ADIC Ignore */
548416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
549416022c9SBarry Smith {
550416022c9SBarry Smith   Mat         A = (Mat) obj;
551416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
552bcd2baecSBarry Smith   ViewerType  vtype;
553bcd2baecSBarry Smith   int         ierr;
554416022c9SBarry Smith 
555bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
556bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
557416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
558416022c9SBarry Smith   }
559bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
560416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
561416022c9SBarry Smith   }
562bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
563416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
564416022c9SBarry Smith   }
565bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
566bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
56717ab2063SBarry Smith   }
56817ab2063SBarry Smith   return 0;
56917ab2063SBarry Smith }
57019bcc07fSBarry Smith 
571c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
5725615d1e5SSatish Balay #undef __FUNC__
5735615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ"
574416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
57517ab2063SBarry Smith {
576416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
57741c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
57894a9d846SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax;
579416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
58017ab2063SBarry Smith 
5816d4a8577SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
58217ab2063SBarry Smith 
58394a9d846SBarry Smith   rmax = ailen[0]; /* determine row with most nonzeros */
58417ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
585416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
58617ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
58794a9d846SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
58817ab2063SBarry Smith     if (fshift) {
589416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
59017ab2063SBarry Smith       N = ailen[i];
59117ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
59217ab2063SBarry Smith         ip[j-fshift] = ip[j];
59317ab2063SBarry Smith         ap[j-fshift] = ap[j];
59417ab2063SBarry Smith       }
59517ab2063SBarry Smith     }
59617ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
59717ab2063SBarry Smith   }
59817ab2063SBarry Smith   if (m) {
59917ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
60017ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
60117ab2063SBarry Smith   }
60217ab2063SBarry Smith   /* reset ilen and imax for each row */
60317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
60417ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
60517ab2063SBarry Smith   }
606416022c9SBarry Smith   a->nz = ai[m] + shift;
60717ab2063SBarry Smith 
60817ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
609416022c9SBarry Smith   if (fshift && a->diag) {
6100452661fSBarry Smith     PetscFree(a->diag);
611416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
612416022c9SBarry Smith     a->diag = 0;
61317ab2063SBarry Smith   }
6144e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
6154e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
6164e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
617b810aeb4SBarry Smith            a->reallocs);
61894a9d846SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax);
619dd5f02e7SSatish Balay   a->reallocs          = 0;
6204e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
6214e220ebcSLois Curfman McInnes 
62276dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
62341c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
62417ab2063SBarry Smith   return 0;
62517ab2063SBarry Smith }
62617ab2063SBarry Smith 
6275615d1e5SSatish Balay #undef __FUNC__
6285615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ"
629416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
63017ab2063SBarry Smith {
631416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
632cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
63317ab2063SBarry Smith   return 0;
63417ab2063SBarry Smith }
635416022c9SBarry Smith 
6365615d1e5SSatish Balay #undef __FUNC__
6375615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ"
63817ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
63917ab2063SBarry Smith {
640416022c9SBarry Smith   Mat        A  = (Mat) obj;
641416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
64290f02eecSBarry Smith   int        ierr;
643d5d45c9bSBarry Smith 
64417ab2063SBarry Smith #if defined(PETSC_LOG)
645416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
64617ab2063SBarry Smith #endif
6470452661fSBarry Smith   PetscFree(a->a);
6480452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
6490452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
6500452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
6510452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
6520452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
65376dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
6540452661fSBarry Smith   PetscFree(a);
65590f02eecSBarry Smith   if (A->mapping) {
65690f02eecSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr);
65790f02eecSBarry Smith   }
658f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
659f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
66017ab2063SBarry Smith   return 0;
66117ab2063SBarry Smith }
66217ab2063SBarry Smith 
6635615d1e5SSatish Balay #undef __FUNC__
6645615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
665416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
66617ab2063SBarry Smith {
66717ab2063SBarry Smith   return 0;
66817ab2063SBarry Smith }
66917ab2063SBarry Smith 
6705615d1e5SSatish Balay #undef __FUNC__
6715615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
672416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
67317ab2063SBarry Smith {
674416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
6756d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)               a->roworiented = 1;
6766d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)            a->roworiented = 0;
6776d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)             a->sorted      = 1;
678219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)           a->sorted      = 0;
6796d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)   a->nonew       = 1;
680*c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew       = -1;
6816d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)  a->nonew       = 0;
6826d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
683219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
6846d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
6856d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
68690f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
68790f02eecSBarry Smith            op == MAT_IGNORE_OFF_PROCESSOR_ENTRIES)
68894a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
6896d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
690e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");}
6916d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
6926d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
6936d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
6946d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
6956d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
696e2f28af5SBarry Smith   else
697e3372554SBarry Smith     {SETERRQ(PETSC_ERR_SUP,0,"unknown option");}
69817ab2063SBarry Smith   return 0;
69917ab2063SBarry Smith }
70017ab2063SBarry Smith 
7015615d1e5SSatish Balay #undef __FUNC__
7025615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
703416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
70417ab2063SBarry Smith {
705416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
706416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
70717ab2063SBarry Smith   Scalar     *x, zero = 0.0;
70817ab2063SBarry Smith 
70917ab2063SBarry Smith   VecSet(&zero,v);
71090f02eecSBarry Smith   VecGetArray_Fast(v,x); VecGetLocalSize(v,&n);
711e3372554SBarry Smith   if (n != a->m) SETERRQ(1,0,"Nonconforming matrix and vector");
712416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
713416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
714416022c9SBarry Smith       if (a->j[j]+shift == i) {
715416022c9SBarry Smith         x[i] = a->a[j];
71617ab2063SBarry Smith         break;
71717ab2063SBarry Smith       }
71817ab2063SBarry Smith     }
71917ab2063SBarry Smith   }
72017ab2063SBarry Smith   return 0;
72117ab2063SBarry Smith }
72217ab2063SBarry Smith 
72317ab2063SBarry Smith /* -------------------------------------------------------*/
72417ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
72517ab2063SBarry Smith /* -------------------------------------------------------*/
7265615d1e5SSatish Balay #undef __FUNC__
7275615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
72844cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
72917ab2063SBarry Smith {
730416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
73117ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
732416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
73317ab2063SBarry Smith 
73490f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
735cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
73617ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
73717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
738416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
739416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
740416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
74117ab2063SBarry Smith     alpha = x[i];
74217ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
74317ab2063SBarry Smith   }
744416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
74517ab2063SBarry Smith   return 0;
74617ab2063SBarry Smith }
747d5d45c9bSBarry Smith 
7485615d1e5SSatish Balay #undef __FUNC__
7495615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
75044cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
75117ab2063SBarry Smith {
752416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
75317ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
754416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
75517ab2063SBarry Smith 
75690f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
75717ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
75817ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
75917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
760416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
761416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
762416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
76317ab2063SBarry Smith     alpha = x[i];
76417ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
76517ab2063SBarry Smith   }
76690f02eecSBarry Smith   PLogFlops(2*a->nz);
76717ab2063SBarry Smith   return 0;
76817ab2063SBarry Smith }
76917ab2063SBarry Smith 
7705615d1e5SSatish Balay #undef __FUNC__
7715615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
77244cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
77317ab2063SBarry Smith {
774416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
77517ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
776416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
77717ab2063SBarry Smith 
77890f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y);
77917ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
780416022c9SBarry Smith   idx  = a->j;
781416022c9SBarry Smith   v    = a->a;
782416022c9SBarry Smith   ii   = a->i;
78317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
784416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
78517ab2063SBarry Smith     sum  = 0.0;
78617ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
78717ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
788416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
78917ab2063SBarry Smith     y[i] = sum;
79017ab2063SBarry Smith   }
791416022c9SBarry Smith   PLogFlops(2*a->nz - m);
79217ab2063SBarry Smith   return 0;
79317ab2063SBarry Smith }
79417ab2063SBarry Smith 
7955615d1e5SSatish Balay #undef __FUNC__
7965615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
79744cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
79817ab2063SBarry Smith {
799416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
80017ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
801cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
80217ab2063SBarry Smith 
80390f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); VecGetArray_Fast(zz,z);
80417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
805cddf8d76SBarry Smith   idx  = a->j;
806cddf8d76SBarry Smith   v    = a->a;
807cddf8d76SBarry Smith   ii   = a->i;
80817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
809cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
81017ab2063SBarry Smith     sum  = y[i];
811cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
812cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
81317ab2063SBarry Smith     z[i] = sum;
81417ab2063SBarry Smith   }
815416022c9SBarry Smith   PLogFlops(2*a->nz);
81617ab2063SBarry Smith   return 0;
81717ab2063SBarry Smith }
81817ab2063SBarry Smith 
81917ab2063SBarry Smith /*
82017ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
82117ab2063SBarry Smith */
82217ab2063SBarry Smith 
8235615d1e5SSatish Balay #undef __FUNC__
8245615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
825416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
82617ab2063SBarry Smith {
827416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
828416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
82917ab2063SBarry Smith 
8300452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
831416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
832416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
833416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
834416022c9SBarry Smith       if (a->j[j]+shift == i) {
83517ab2063SBarry Smith         diag[i] = j - shift;
83617ab2063SBarry Smith         break;
83717ab2063SBarry Smith       }
83817ab2063SBarry Smith     }
83917ab2063SBarry Smith   }
840416022c9SBarry Smith   a->diag = diag;
84117ab2063SBarry Smith   return 0;
84217ab2063SBarry Smith }
84317ab2063SBarry Smith 
8445615d1e5SSatish Balay #undef __FUNC__
8455615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
84644cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
84717ab2063SBarry Smith                            double fshift,int its,Vec xx)
84817ab2063SBarry Smith {
849416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
850416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
851d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
85217ab2063SBarry Smith 
85390f02eecSBarry Smith   VecGetArray_Fast(xx,x); VecGetArray_Fast(bb,b);
854416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
855416022c9SBarry Smith   diag = a->diag;
856416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
85717ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
85817ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
85917ab2063SBarry Smith     bs = b + shift;
86017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
861416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
862416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
863416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
864416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
86517ab2063SBarry Smith         sum  = b[i]*d/omega;
86617ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
86717ab2063SBarry Smith         x[i] = sum;
86817ab2063SBarry Smith     }
86917ab2063SBarry Smith     return 0;
87017ab2063SBarry Smith   }
87117ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
872e3372554SBarry Smith     SETERRQ(1,0,"SOR_APPLY_LOWER is not done");
87317ab2063SBarry Smith   }
874416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
87517ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
87617ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
87717ab2063SBarry Smith 
87817ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
87917ab2063SBarry Smith 
88017ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
88117ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
88217ab2063SBarry Smith     is the relaxation factor.
88317ab2063SBarry Smith     */
8840452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
88517ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
88617ab2063SBarry Smith 
88717ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
88817ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
889416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
890416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
891416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
892416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
89317ab2063SBarry Smith       sum  = b[i];
89417ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
89517ab2063SBarry Smith       x[i] = omega*(sum/d);
89617ab2063SBarry Smith     }
89717ab2063SBarry Smith 
89817ab2063SBarry Smith     /*  t = b - (2*E - D)x */
899416022c9SBarry Smith     v = a->a;
90017ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
90117ab2063SBarry Smith 
90217ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
903416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
904416022c9SBarry Smith     diag = a->diag;
90517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
906416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
907416022c9SBarry Smith       n    = diag[i] - a->i[i];
908416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
909416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
91017ab2063SBarry Smith       sum  = t[i];
91117ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
91217ab2063SBarry Smith       t[i] = omega*(sum/d);
91317ab2063SBarry Smith     }
91417ab2063SBarry Smith 
91517ab2063SBarry Smith     /*  x = x + t */
91617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
9170452661fSBarry Smith     PetscFree(t);
91817ab2063SBarry Smith     return 0;
91917ab2063SBarry Smith   }
92017ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
92117ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
92217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
923416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
924416022c9SBarry Smith         n    = diag[i] - a->i[i];
925416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
926416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
92717ab2063SBarry Smith         sum  = b[i];
92817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
92917ab2063SBarry Smith         x[i] = omega*(sum/d);
93017ab2063SBarry Smith       }
93117ab2063SBarry Smith       xb = x;
93217ab2063SBarry Smith     }
93317ab2063SBarry Smith     else xb = b;
93417ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
93517ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
93617ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
937416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
93817ab2063SBarry Smith       }
93917ab2063SBarry Smith     }
94017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
94117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
942416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
943416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
944416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
945416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
94617ab2063SBarry Smith         sum  = xb[i];
94717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
94817ab2063SBarry Smith         x[i] = omega*(sum/d);
94917ab2063SBarry Smith       }
95017ab2063SBarry Smith     }
95117ab2063SBarry Smith     its--;
95217ab2063SBarry Smith   }
95317ab2063SBarry Smith   while (its--) {
95417ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
95517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
956416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
957416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
958416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
959416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
96017ab2063SBarry Smith         sum  = b[i];
96117ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
9627e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
96317ab2063SBarry Smith       }
96417ab2063SBarry Smith     }
96517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
96617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
967416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
968416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
969416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
970416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
97117ab2063SBarry Smith         sum  = b[i];
97217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
9737e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
97417ab2063SBarry Smith       }
97517ab2063SBarry Smith     }
97617ab2063SBarry Smith   }
97717ab2063SBarry Smith   return 0;
97817ab2063SBarry Smith }
97917ab2063SBarry Smith 
9805615d1e5SSatish Balay #undef __FUNC__
9815615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
9824e220ebcSLois Curfman McInnes static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
98317ab2063SBarry Smith {
984416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
9854e220ebcSLois Curfman McInnes 
9864e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
9874e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
9884e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
9894e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
9904e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
9914e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
9924e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
9934e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
9944e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
9954e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
9964e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
9974e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
9984e220ebcSLois Curfman McInnes   info->memory         = A->mem;
9994e220ebcSLois Curfman McInnes   if (A->factor) {
10004e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10014e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
10024e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
10034e220ebcSLois Curfman McInnes   } else {
10044e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
10054e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
10064e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
10074e220ebcSLois Curfman McInnes   }
100817ab2063SBarry Smith   return 0;
100917ab2063SBarry Smith }
101017ab2063SBarry Smith 
101117ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
101217ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
101317ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
101417ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
101517ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
101617ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
101717ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
101817ab2063SBarry Smith 
10195615d1e5SSatish Balay #undef __FUNC__
10205615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
102117ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
102217ab2063SBarry Smith {
1023416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1024416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
102517ab2063SBarry Smith 
102677c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
102717ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
102817ab2063SBarry Smith   if (diag) {
102917ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1030e3372554SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range");
1031416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1032416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1033416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1034416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
103517ab2063SBarry Smith       }
103617ab2063SBarry Smith       else {
103717ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
103817ab2063SBarry Smith         CHKERRQ(ierr);
103917ab2063SBarry Smith       }
104017ab2063SBarry Smith     }
104117ab2063SBarry Smith   }
104217ab2063SBarry Smith   else {
104317ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1044e3372554SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range");
1045416022c9SBarry Smith       a->ilen[rows[i]] = 0;
104617ab2063SBarry Smith     }
104717ab2063SBarry Smith   }
104817ab2063SBarry Smith   ISRestoreIndices(is,&rows);
104943a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
105017ab2063SBarry Smith   return 0;
105117ab2063SBarry Smith }
105217ab2063SBarry Smith 
10535615d1e5SSatish Balay #undef __FUNC__
10545615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
1055416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
105617ab2063SBarry Smith {
1057416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1058416022c9SBarry Smith   *m = a->m; *n = a->n;
105917ab2063SBarry Smith   return 0;
106017ab2063SBarry Smith }
106117ab2063SBarry Smith 
10625615d1e5SSatish Balay #undef __FUNC__
10635615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
1064416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
106517ab2063SBarry Smith {
1066416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1067416022c9SBarry Smith   *m = 0; *n = a->m;
106817ab2063SBarry Smith   return 0;
106917ab2063SBarry Smith }
1070026e39d0SSatish Balay 
10715615d1e5SSatish Balay #undef __FUNC__
10725615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
10734e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
107417ab2063SBarry Smith {
1075416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1076c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
107717ab2063SBarry Smith 
1078e3372554SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,0,"Row out of range");
107917ab2063SBarry Smith 
1080416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1081416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
108217ab2063SBarry Smith   if (idx) {
1083416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
10844e093b46SBarry Smith     if (*nz && shift) {
10850452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
108617ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
10874e093b46SBarry Smith     } else if (*nz) {
10884e093b46SBarry Smith       *idx = itmp;
108917ab2063SBarry Smith     }
109017ab2063SBarry Smith     else *idx = 0;
109117ab2063SBarry Smith   }
109217ab2063SBarry Smith   return 0;
109317ab2063SBarry Smith }
109417ab2063SBarry Smith 
10955615d1e5SSatish Balay #undef __FUNC__
10965615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
10974e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
109817ab2063SBarry Smith {
10994e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11004e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
110117ab2063SBarry Smith   return 0;
110217ab2063SBarry Smith }
110317ab2063SBarry Smith 
11045615d1e5SSatish Balay #undef __FUNC__
11055615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
1106cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
110717ab2063SBarry Smith {
1108416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1109416022c9SBarry Smith   Scalar     *v = a->a;
111017ab2063SBarry Smith   double     sum = 0.0;
1111416022c9SBarry Smith   int        i, j,shift = a->indexshift;
111217ab2063SBarry Smith 
111317ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1114416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
111517ab2063SBarry Smith #if defined(PETSC_COMPLEX)
111617ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
111717ab2063SBarry Smith #else
111817ab2063SBarry Smith       sum += (*v)*(*v); v++;
111917ab2063SBarry Smith #endif
112017ab2063SBarry Smith     }
112117ab2063SBarry Smith     *norm = sqrt(sum);
112217ab2063SBarry Smith   }
112317ab2063SBarry Smith   else if (type == NORM_1) {
112417ab2063SBarry Smith     double *tmp;
1125416022c9SBarry Smith     int    *jj = a->j;
11260452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
1127cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
112817ab2063SBarry Smith     *norm = 0.0;
1129416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1130a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
113117ab2063SBarry Smith     }
1132416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
113317ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
113417ab2063SBarry Smith     }
11350452661fSBarry Smith     PetscFree(tmp);
113617ab2063SBarry Smith   }
113717ab2063SBarry Smith   else if (type == NORM_INFINITY) {
113817ab2063SBarry Smith     *norm = 0.0;
1139416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1140416022c9SBarry Smith       v = a->a + a->i[j] + shift;
114117ab2063SBarry Smith       sum = 0.0;
1142416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1143cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
114417ab2063SBarry Smith       }
114517ab2063SBarry Smith       if (sum > *norm) *norm = sum;
114617ab2063SBarry Smith     }
114717ab2063SBarry Smith   }
114817ab2063SBarry Smith   else {
1149e3372554SBarry Smith     SETERRQ(1,0,"No support for two norm yet");
115017ab2063SBarry Smith   }
115117ab2063SBarry Smith   return 0;
115217ab2063SBarry Smith }
115317ab2063SBarry Smith 
11545615d1e5SSatish Balay #undef __FUNC__
11555615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
1156416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
115717ab2063SBarry Smith {
1158416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1159416022c9SBarry Smith   Mat        C;
1160416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1161416022c9SBarry Smith   int        shift = a->indexshift;
1162d5d45c9bSBarry Smith   Scalar     *array = a->a;
116317ab2063SBarry Smith 
11643638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
1165e3372554SBarry Smith     SETERRQ(1,0,"Square matrix only for in-place");
11660452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1167cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
116817ab2063SBarry Smith   if (shift) {
116917ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
117017ab2063SBarry Smith   }
117117ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1172416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
11730452661fSBarry Smith   PetscFree(col);
117417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
117517ab2063SBarry Smith     len = ai[i+1]-ai[i];
1176416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
117717ab2063SBarry Smith     array += len; aj += len;
117817ab2063SBarry Smith   }
117917ab2063SBarry Smith   if (shift) {
118017ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
118117ab2063SBarry Smith   }
118217ab2063SBarry Smith 
11836d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11846d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
118517ab2063SBarry Smith 
11863638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1187416022c9SBarry Smith     *B = C;
118817ab2063SBarry Smith   } else {
1189416022c9SBarry Smith     /* This isn't really an in-place transpose */
11900452661fSBarry Smith     PetscFree(a->a);
11910452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
11920452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
11930452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
11940452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
11950452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
11961073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
11970452661fSBarry Smith     PetscFree(a);
1198416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
11990452661fSBarry Smith     PetscHeaderDestroy(C);
120017ab2063SBarry Smith   }
120117ab2063SBarry Smith   return 0;
120217ab2063SBarry Smith }
120317ab2063SBarry Smith 
12045615d1e5SSatish Balay #undef __FUNC__
12055615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
1206f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
120717ab2063SBarry Smith {
1208416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
120917ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1210d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
121117ab2063SBarry Smith 
121217ab2063SBarry Smith   if (ll) {
12133ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
12143ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
12159b1297e1SSatish Balay     VecGetLocalSize_Fast(ll,m);
1216e3372554SBarry Smith     if (m != a->m) SETERRQ(1,0,"Left scaling vector wrong length");
121790f02eecSBarry Smith     VecGetArray_Fast(ll,l);
1218416022c9SBarry Smith     v = a->a;
121917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
122017ab2063SBarry Smith       x = l[i];
1221416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
122217ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
122317ab2063SBarry Smith     }
122444cd7ae7SLois Curfman McInnes     PLogFlops(nz);
122517ab2063SBarry Smith   }
122617ab2063SBarry Smith   if (rr) {
12279b1297e1SSatish Balay     VecGetLocalSize_Fast(rr,n);
1228e3372554SBarry Smith     if (n != a->n) SETERRQ(1,0,"Right scaling vector wrong length");
122990f02eecSBarry Smith     VecGetArray_Fast(rr,r);
1230416022c9SBarry Smith     v = a->a; jj = a->j;
123117ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
123217ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
123317ab2063SBarry Smith     }
123444cd7ae7SLois Curfman McInnes     PLogFlops(nz);
123517ab2063SBarry Smith   }
123617ab2063SBarry Smith   return 0;
123717ab2063SBarry Smith }
123817ab2063SBarry Smith 
12395615d1e5SSatish Balay #undef __FUNC__
12405615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
1241cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
124217ab2063SBarry Smith {
1243db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
124402834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
124599141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1246a2744918SBarry Smith   register int sum,lensi;
124799141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
124899141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
124999141d43SSatish Balay   Scalar       *a_new,*mat_a;
1250416022c9SBarry Smith   Mat          C;
125117ab2063SBarry Smith 
1252b48a1e75SSatish Balay   ierr = ISSorted(isrow,(PetscTruth*)&i);
1253e3372554SBarry Smith   if (!i) SETERRQ(1,0,"ISrow is not sorted");
125499141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
1255e3372554SBarry Smith   if (!i) SETERRQ(1,0,"IScol is not sorted");
125699141d43SSatish Balay 
125717ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
125817ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
125917ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
126017ab2063SBarry Smith 
12617264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
126202834360SBarry Smith     /* special case of contiguous rows */
126357faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
126402834360SBarry Smith     starts = lens + ncols;
126502834360SBarry Smith     /* loop over new rows determining lens and starting points */
126602834360SBarry Smith     for (i=0; i<nrows; i++) {
1267a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1268a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
126902834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1270d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
127102834360SBarry Smith           starts[i] = k;
127202834360SBarry Smith           break;
127302834360SBarry Smith 	}
127402834360SBarry Smith       }
1275a2744918SBarry Smith       sum = 0;
127602834360SBarry Smith       while (k < kend) {
1277d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1278a2744918SBarry Smith         sum++;
127902834360SBarry Smith       }
1280a2744918SBarry Smith       lens[i] = sum;
128102834360SBarry Smith     }
128202834360SBarry Smith     /* create submatrix */
1283cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
128408480c60SBarry Smith       int n_cols,n_rows;
128508480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1286e3372554SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,0,"");
1287d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
128808480c60SBarry Smith       C = *B;
128908480c60SBarry Smith     }
129008480c60SBarry Smith     else {
129102834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
129208480c60SBarry Smith     }
1293db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1294db02288aSLois Curfman McInnes 
129502834360SBarry Smith     /* loop over rows inserting into submatrix */
1296db02288aSLois Curfman McInnes     a_new    = c->a;
1297db02288aSLois Curfman McInnes     j_new    = c->j;
1298db02288aSLois Curfman McInnes     i_new    = c->i;
1299db02288aSLois Curfman McInnes     i_new[0] = -shift;
130002834360SBarry Smith     for (i=0; i<nrows; i++) {
1301a2744918SBarry Smith       ii    = starts[i];
1302a2744918SBarry Smith       lensi = lens[i];
1303a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1304a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
130502834360SBarry Smith       }
1306a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1307a2744918SBarry Smith       a_new      += lensi;
1308a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1309a2744918SBarry Smith       c->ilen[i]  = lensi;
131002834360SBarry Smith     }
13110452661fSBarry Smith     PetscFree(lens);
131202834360SBarry Smith   }
131302834360SBarry Smith   else {
131402834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
13150452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
131602834360SBarry Smith     ssmap = smap + shift;
131799141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1318cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
131917ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
132002834360SBarry Smith     /* determine lens of each row */
132102834360SBarry Smith     for (i=0; i<nrows; i++) {
1322d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
132302834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
132402834360SBarry Smith       lens[i] = 0;
132502834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1326d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
132702834360SBarry Smith           lens[i]++;
132802834360SBarry Smith         }
132902834360SBarry Smith       }
133002834360SBarry Smith     }
133117ab2063SBarry Smith     /* Create and fill new matrix */
1332a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
133399141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
133499141d43SSatish Balay 
1335e3372554SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(1,0,"");
133699141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1337e3372554SBarry Smith         SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros");
133899141d43SSatish Balay       }
133999141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
134008480c60SBarry Smith       C = *B;
134108480c60SBarry Smith     }
134208480c60SBarry Smith     else {
134302834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
134408480c60SBarry Smith     }
134599141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
134617ab2063SBarry Smith     for (i=0; i<nrows; i++) {
134799141d43SSatish Balay       row    = irow[i];
134817ab2063SBarry Smith       nznew  = 0;
134999141d43SSatish Balay       kstart = ai[row]+shift;
135099141d43SSatish Balay       kend   = kstart + a->ilen[row];
135199141d43SSatish Balay       mat_i  = c->i[i]+shift;
135299141d43SSatish Balay       mat_j  = c->j + mat_i;
135399141d43SSatish Balay       mat_a  = c->a + mat_i;
135499141d43SSatish Balay       mat_ilen = c->ilen + i;
135517ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
135699141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
135799141d43SSatish Balay           *mat_j++ = tcol - (!shift);
135899141d43SSatish Balay           *mat_a++ = a->a[k];
135999141d43SSatish Balay           (*mat_ilen)++;
136099141d43SSatish Balay 
136117ab2063SBarry Smith         }
136217ab2063SBarry Smith       }
136317ab2063SBarry Smith     }
136402834360SBarry Smith     /* Free work space */
136502834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
136699141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
136702834360SBarry Smith   }
13686d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13696d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
137017ab2063SBarry Smith 
137117ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1372416022c9SBarry Smith   *B = C;
137317ab2063SBarry Smith   return 0;
137417ab2063SBarry Smith }
137517ab2063SBarry Smith 
1376a871dcd8SBarry Smith /*
137763b91edcSBarry Smith      note: This can only work for identity for row and col. It would
137863b91edcSBarry Smith    be good to check this and otherwise generate an error.
1379a871dcd8SBarry Smith */
13805615d1e5SSatish Balay #undef __FUNC__
13815615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
138263b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1383a871dcd8SBarry Smith {
138463b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
138508480c60SBarry Smith   int        ierr;
138663b91edcSBarry Smith   Mat        outA;
138763b91edcSBarry Smith 
1388e3372554SBarry Smith   if (fill != 0) SETERRQ(1,0,"Only fill=0 supported");
1389a871dcd8SBarry Smith 
139063b91edcSBarry Smith   outA          = inA;
139163b91edcSBarry Smith   inA->factor   = FACTOR_LU;
139263b91edcSBarry Smith   a->row        = row;
139363b91edcSBarry Smith   a->col        = col;
139463b91edcSBarry Smith 
139594a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
13960452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
139794a9d846SBarry Smith   }
139863b91edcSBarry Smith 
139908480c60SBarry Smith   if (!a->diag) {
140008480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
140163b91edcSBarry Smith   }
140263b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1403a871dcd8SBarry Smith   return 0;
1404a871dcd8SBarry Smith }
1405a871dcd8SBarry Smith 
1406f0b747eeSBarry Smith #include "pinclude/plapack.h"
14075615d1e5SSatish Balay #undef __FUNC__
14085615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
1409f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1410f0b747eeSBarry Smith {
1411f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1412f0b747eeSBarry Smith   int        one = 1;
1413f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1414f0b747eeSBarry Smith   PLogFlops(a->nz);
1415f0b747eeSBarry Smith   return 0;
1416f0b747eeSBarry Smith }
1417f0b747eeSBarry Smith 
14185615d1e5SSatish Balay #undef __FUNC__
14195615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
1420cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1421cddf8d76SBarry Smith                                     Mat **B)
1422cddf8d76SBarry Smith {
1423cddf8d76SBarry Smith   int ierr,i;
1424cddf8d76SBarry Smith 
1425cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
14260452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1427cddf8d76SBarry Smith   }
1428cddf8d76SBarry Smith 
1429cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1430905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1431cddf8d76SBarry Smith   }
1432cddf8d76SBarry Smith   return 0;
1433cddf8d76SBarry Smith }
1434cddf8d76SBarry Smith 
14355615d1e5SSatish Balay #undef __FUNC__
14365615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
14375a838052SSatish Balay static int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
14385a838052SSatish Balay {
14395a838052SSatish Balay   *bs = 1;
14405a838052SSatish Balay   return 0;
14415a838052SSatish Balay }
14425a838052SSatish Balay 
14435615d1e5SSatish Balay #undef __FUNC__
14445615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
1445e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
14464dcbc457SBarry Smith {
1447e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
144806763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
14498a047759SSatish Balay   int        start, end, *ai, *aj;
145006763907SSatish Balay   char       *table;
14518a047759SSatish Balay   shift = a->indexshift;
1452e4d965acSSatish Balay   m     = a->m;
1453e4d965acSSatish Balay   ai    = a->i;
14548a047759SSatish Balay   aj    = a->j+shift;
14558a047759SSatish Balay 
1456e3372554SBarry Smith   if (ov < 0)  SETERRQ(1,0,"illegal overlap value used");
145706763907SSatish Balay 
145806763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
145906763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
146006763907SSatish Balay 
1461e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1462b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1463e4d965acSSatish Balay     isz  = 0;
146406763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1465e4d965acSSatish Balay 
1466e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
14674dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
146877c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1469e4d965acSSatish Balay 
1470dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1471e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
147206763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
14734dcbc457SBarry Smith     }
147406763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
147506763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1476e4d965acSSatish Balay 
147704a348a9SBarry Smith     k = 0;
147804a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
147904a348a9SBarry Smith       n = isz;
148006763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1481e4d965acSSatish Balay         row   = nidx[k];
1482e4d965acSSatish Balay         start = ai[row];
1483e4d965acSSatish Balay         end   = ai[row+1];
148404a348a9SBarry Smith         for ( l = start; l<end ; l++){
14858a047759SSatish Balay           val = aj[l] + shift;
148606763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1487e4d965acSSatish Balay         }
1488e4d965acSSatish Balay       }
1489e4d965acSSatish Balay     }
1490537820f0SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1491e4d965acSSatish Balay   }
149204a348a9SBarry Smith   PetscFree(table);
149306763907SSatish Balay   PetscFree(nidx);
1494e4d965acSSatish Balay   return 0;
14954dcbc457SBarry Smith }
149617ab2063SBarry Smith 
14970513a670SBarry Smith /* -------------------------------------------------------------- */
14980513a670SBarry Smith #undef __FUNC__
14990513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
15000513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
15010513a670SBarry Smith {
15020513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
15030513a670SBarry Smith   Scalar     *vwork;
15040513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
15050513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
15060513a670SBarry Smith 
15070513a670SBarry Smith   ierr = ISGetIndices(rowp,&row); CHKERRQ(ierr);
15080513a670SBarry Smith   ierr = ISGetIndices(colp,&col); CHKERRQ(ierr);
15090513a670SBarry Smith 
15100513a670SBarry Smith   /* determine lengths of permuted rows */
15110513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
15120513a670SBarry Smith   for (i=0; i<m; i++ ) {
15130513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
15140513a670SBarry Smith   }
15150513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
15160513a670SBarry Smith   PetscFree(lens);
15170513a670SBarry Smith 
15180513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
15190513a670SBarry Smith   for (i=0; i<m; i++) {
15200513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15210513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
15220513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
15230513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
15240513a670SBarry Smith   }
15250513a670SBarry Smith   PetscFree(cnew);
15260513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
15270513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
15280513a670SBarry Smith   ierr = ISRestoreIndices(rowp,&row); CHKERRQ(ierr);
15290513a670SBarry Smith   ierr = ISRestoreIndices(colp,&col); CHKERRQ(ierr);
15300513a670SBarry Smith   return 0;
15310513a670SBarry Smith }
15320513a670SBarry Smith 
15335615d1e5SSatish Balay #undef __FUNC__
15345615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1535682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1536682d7d0cSBarry Smith {
1537682d7d0cSBarry Smith   static int called = 0;
1538682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1539682d7d0cSBarry Smith 
1540682d7d0cSBarry Smith   if (called) return 0; else called = 1;
154177c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
15420f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
15430f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
15440f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_no_inode: Do not use inodes\n");
15450f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1546682d7d0cSBarry Smith #if defined(HAVE_ESSL)
15470f665d81SLois Curfman McInnes   PetscPrintf(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1548682d7d0cSBarry Smith #endif
1549682d7d0cSBarry Smith   return 0;
1550682d7d0cSBarry Smith }
1551205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1552a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1553a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1554a93ec695SBarry Smith 
1555682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
155617ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
155717ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1558416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1559416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
156017ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
156117ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
156217ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
156317ab2063SBarry Smith        MatRelax_SeqAIJ,
156417ab2063SBarry Smith        MatTranspose_SeqAIJ,
15657264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1566f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
156717ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
156817ab2063SBarry Smith        MatCompress_SeqAIJ,
156917ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
157017ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
157117ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
157217ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
157394a9d846SBarry Smith        0,0,
15743d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1575cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
15767eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1577682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1578f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
15795a838052SSatish Balay        MatScale_SeqAIJ,0,0,
15806945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
15816945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
15823b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
15833b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
15843b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1585a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1586a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
15870513a670SBarry Smith        MatColoringPatch_SeqAIJ,
15880513a670SBarry Smith        0,
15890513a670SBarry Smith        MatPermute_SeqAIJ};
159017ab2063SBarry Smith 
159117ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
159217ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
159317ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
159417ab2063SBarry Smith 
15955615d1e5SSatish Balay #undef __FUNC__
15965615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
159717ab2063SBarry Smith /*@C
1598682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
15990d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
16006e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
16012bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
16022bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
160317ab2063SBarry Smith 
160417ab2063SBarry Smith    Input Parameters:
160517ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
160617ab2063SBarry Smith .  m - number of rows
160717ab2063SBarry Smith .  n - number of columns
160817ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
16092bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
16102bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
161117ab2063SBarry Smith 
161217ab2063SBarry Smith    Output Parameter:
1613416022c9SBarry Smith .  A - the matrix
161417ab2063SBarry Smith 
161517ab2063SBarry Smith    Notes:
161617ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
161717ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
16180002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
161944cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
162017ab2063SBarry Smith 
162117ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1622a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
16233d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
16246da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
162517ab2063SBarry Smith 
1626682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1627682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1628682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
16296c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
16306c7ebb05SLois Curfman McInnes 
16316c7ebb05SLois Curfman McInnes    Options Database Keys:
16326c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
16336e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
16346e62573dSLois Curfman McInnes $        (max limit=5)
16356e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
16366e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
16376e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
163817ab2063SBarry Smith 
163917ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
164017ab2063SBarry Smith @*/
1641416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
164217ab2063SBarry Smith {
1643416022c9SBarry Smith   Mat        B;
1644416022c9SBarry Smith   Mat_SeqAIJ *b;
16456945ee14SBarry Smith   int        i, len, ierr, flg,size;
16466945ee14SBarry Smith 
16476945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1648e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"Comm must be of size 1");
1649d5d45c9bSBarry Smith 
1650416022c9SBarry Smith   *A                  = 0;
16510452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1652416022c9SBarry Smith   PLogObjectCreate(B);
16530452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
165444cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1655416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1656416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1657416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1658416022c9SBarry Smith   B->factor           = 0;
1659416022c9SBarry Smith   B->lupivotthreshold = 1.0;
166090f02eecSBarry Smith   B->mapping          = 0;
16617a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
166269957df2SSatish Balay                           &flg); CHKERRQ(ierr);
16637a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
16647a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
16657a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1666416022c9SBarry Smith   b->row              = 0;
1667416022c9SBarry Smith   b->col              = 0;
1668416022c9SBarry Smith   b->indexshift       = 0;
1669b810aeb4SBarry Smith   b->reallocs         = 0;
167069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
167169957df2SSatish Balay   if (flg) b->indexshift = -1;
167217ab2063SBarry Smith 
167344cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
167444cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
16750452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1676b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
16777b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
16787b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1679416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
168017ab2063SBarry Smith     nz = nz*m;
168117ab2063SBarry Smith   }
168217ab2063SBarry Smith   else {
168317ab2063SBarry Smith     nz = 0;
1684416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
168517ab2063SBarry Smith   }
168617ab2063SBarry Smith 
168717ab2063SBarry Smith   /* allocate the matrix space */
1688416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
16890452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1690416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1691cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1692416022c9SBarry Smith   b->i  = b->j + nz;
1693416022c9SBarry Smith   b->singlemalloc = 1;
169417ab2063SBarry Smith 
1695416022c9SBarry Smith   b->i[0] = -b->indexshift;
169617ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1697416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
169817ab2063SBarry Smith   }
169917ab2063SBarry Smith 
1700416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
17010452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1702416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1703416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
170417ab2063SBarry Smith 
1705416022c9SBarry Smith   b->nz               = 0;
1706416022c9SBarry Smith   b->maxnz            = nz;
1707416022c9SBarry Smith   b->sorted           = 0;
1708416022c9SBarry Smith   b->roworiented      = 1;
1709416022c9SBarry Smith   b->nonew            = 0;
1710416022c9SBarry Smith   b->diag             = 0;
1711416022c9SBarry Smith   b->solve_work       = 0;
171271bd300dSLois Curfman McInnes   b->spptr            = 0;
1713754ec7b1SSatish Balay   b->inode.node_count = 0;
1714754ec7b1SSatish Balay   b->inode.size       = 0;
17156c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
17166c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
17174e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
171817ab2063SBarry Smith 
1719416022c9SBarry Smith   *A = B;
17204e220ebcSLois Curfman McInnes 
17214b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
17224b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
172369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
172469957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
17254b14c69eSBarry Smith #endif
172669957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
172769957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
172869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
172969957df2SSatish Balay   if (flg) {
1730e3372554SBarry Smith     if (!b->indexshift) SETERRQ(1,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1731416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
173217ab2063SBarry Smith   }
173369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
173469957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
173517ab2063SBarry Smith   return 0;
173617ab2063SBarry Smith }
173717ab2063SBarry Smith 
17385615d1e5SSatish Balay #undef __FUNC__
17395615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
17403d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
174117ab2063SBarry Smith {
1742416022c9SBarry Smith   Mat        C;
1743416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
174408480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
174517ab2063SBarry Smith 
17464043dd9cSLois Curfman McInnes   *B = 0;
17470452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1748416022c9SBarry Smith   PLogObjectCreate(C);
17490452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
175041c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1751416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1752416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1753416022c9SBarry Smith   C->factor     = A->factor;
1754416022c9SBarry Smith   c->row        = 0;
1755416022c9SBarry Smith   c->col        = 0;
1756416022c9SBarry Smith   c->indexshift = shift;
1757c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
175817ab2063SBarry Smith 
175944cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
176044cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
176144cd7ae7SLois Curfman McInnes   C->M          = a->m;
176244cd7ae7SLois Curfman McInnes   C->N          = a->n;
176317ab2063SBarry Smith 
17640452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
17650452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
176617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1767416022c9SBarry Smith     c->imax[i] = a->imax[i];
1768416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
176917ab2063SBarry Smith   }
177017ab2063SBarry Smith 
177117ab2063SBarry Smith   /* allocate the matrix space */
1772416022c9SBarry Smith   c->singlemalloc = 1;
1773416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
17740452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1775416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1776416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1777416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
177817ab2063SBarry Smith   if (m > 0) {
1779416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
178008480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1781416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
178217ab2063SBarry Smith     }
178308480c60SBarry Smith   }
178417ab2063SBarry Smith 
1785416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1786416022c9SBarry Smith   c->sorted      = a->sorted;
1787416022c9SBarry Smith   c->roworiented = a->roworiented;
1788416022c9SBarry Smith   c->nonew       = a->nonew;
17897a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
179017ab2063SBarry Smith 
1791416022c9SBarry Smith   if (a->diag) {
17920452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1793416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
179417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1795416022c9SBarry Smith       c->diag[i] = a->diag[i];
179617ab2063SBarry Smith     }
179717ab2063SBarry Smith   }
1798416022c9SBarry Smith   else c->diag          = 0;
17996c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
18006c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1801754ec7b1SSatish Balay   if (a->inode.size){
1802daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1803754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1804daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1805754ec7b1SSatish Balay   } else {
1806754ec7b1SSatish Balay     c->inode.size       = 0;
1807754ec7b1SSatish Balay     c->inode.node_count = 0;
1808754ec7b1SSatish Balay   }
1809416022c9SBarry Smith   c->nz                 = a->nz;
1810416022c9SBarry Smith   c->maxnz              = a->maxnz;
1811416022c9SBarry Smith   c->solve_work         = 0;
181276dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1813754ec7b1SSatish Balay 
1814416022c9SBarry Smith   *B = C;
181517ab2063SBarry Smith   return 0;
181617ab2063SBarry Smith }
181717ab2063SBarry Smith 
18185615d1e5SSatish Balay #undef __FUNC__
18195615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
182019bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
182117ab2063SBarry Smith {
1822416022c9SBarry Smith   Mat_SeqAIJ   *a;
1823416022c9SBarry Smith   Mat          B;
182417699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1825bcd2baecSBarry Smith   MPI_Comm     comm;
182617ab2063SBarry Smith 
182719bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
182817699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
1829e3372554SBarry Smith   if (size > 1) SETERRQ(1,0,"view must have one processor");
183090ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
183177c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1832e3372554SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object in file");
183317ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
183417ab2063SBarry Smith 
183517ab2063SBarry Smith   /* read in row lengths */
18360452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
183777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
183817ab2063SBarry Smith 
183917ab2063SBarry Smith   /* create our matrix */
1840416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1841416022c9SBarry Smith   B = *A;
1842416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1843416022c9SBarry Smith   shift = a->indexshift;
184417ab2063SBarry Smith 
184517ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
184677c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
184717ab2063SBarry Smith   if (shift) {
184817ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1849416022c9SBarry Smith       a->j[i] += 1;
185017ab2063SBarry Smith     }
185117ab2063SBarry Smith   }
185217ab2063SBarry Smith 
185317ab2063SBarry Smith   /* read in nonzero values */
185477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
185517ab2063SBarry Smith 
185617ab2063SBarry Smith   /* set matrix "i" values */
1857416022c9SBarry Smith   a->i[0] = -shift;
185817ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1859416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1860416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
186117ab2063SBarry Smith   }
18620452661fSBarry Smith   PetscFree(rowlengths);
186317ab2063SBarry Smith 
18646d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
18656d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
186617ab2063SBarry Smith   return 0;
186717ab2063SBarry Smith }
186817ab2063SBarry Smith 
18695615d1e5SSatish Balay #undef __FUNC__
18705615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
1871686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
18727264ac53SSatish Balay {
18737264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
18747264ac53SSatish Balay 
1875e3372554SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,0,"Matrices must be same type");
18767264ac53SSatish Balay 
18777264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
18787264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1879bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
188077c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1881bcd2baecSBarry Smith   }
18827264ac53SSatish Balay 
18837264ac53SSatish Balay   /* if the a->i are the same */
1884bcd2baecSBarry Smith   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
188577c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
18867264ac53SSatish Balay   }
18877264ac53SSatish Balay 
18887264ac53SSatish Balay   /* if a->j are the same */
1889bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
189077c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1891bcd2baecSBarry Smith   }
1892bcd2baecSBarry Smith 
1893bcd2baecSBarry Smith   /* if a->a are the same */
189419bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
189577c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
18967264ac53SSatish Balay   }
189777c4ece6SBarry Smith   *flg = PETSC_TRUE;
18987264ac53SSatish Balay   return 0;
18997264ac53SSatish Balay 
19007264ac53SSatish Balay }
1901