xref: /petsc/src/mat/impls/aij/seq/aij.c (revision d8ced48e40de490b1575e35d1e64a4a0df4d5aad)
117ab2063SBarry Smith #ifndef lint
2*d8ced48eSBarry Smith static char vcid[] = "$Id: aij.c,v 1.132 1996/01/04 15:53:20 curfman Exp bsmith $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5d5d45c9bSBarry Smith /*
6d5d45c9bSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7d5d45c9bSBarry Smith   matrix storage format.
8d5d45c9bSBarry Smith */
917ab2063SBarry Smith #include "aij.h"
1017ab2063SBarry Smith #include "vec/vecimpl.h"
1117ab2063SBarry Smith #include "inline/spops.h"
1217ab2063SBarry Smith 
1317ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**);
1417ab2063SBarry Smith 
15416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
1617ab2063SBarry Smith {
17416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
18a2744918SBarry Smith   int        ierr, *ia, *ja,n,*idx,i;
193d54f168SSatish Balay   /*Viewer     V1, V2;*/
2017ab2063SBarry Smith 
21416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix");
2217ab2063SBarry Smith 
23a2744918SBarry Smith   /*
24a2744918SBarry Smith      this is tacky: In the future when we have written special factorization
25a2744918SBarry Smith      and solve routines for the identity permutation we should use a
26a2744918SBarry Smith      stride index set instead of the general one.
27a2744918SBarry Smith   */
28a2744918SBarry Smith   if (type  == ORDER_NATURAL) {
29a2744918SBarry Smith     n = a->n;
30a2744918SBarry Smith     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
31a2744918SBarry Smith     for ( i=0; i<n; i++ ) idx[i] = i;
32a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
33a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
34a2744918SBarry Smith     PetscFree(idx);
35a2744918SBarry Smith     ISSetPermutation(*rperm);
36a2744918SBarry Smith     ISSetPermutation(*cperm);
37a2744918SBarry Smith     ISSetIdentity(*rperm);
38a2744918SBarry Smith     ISSetIdentity(*cperm);
39a2744918SBarry Smith     return 0;
40a2744918SBarry Smith   }
41a2744918SBarry Smith 
42416022c9SBarry Smith   ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr);
43416022c9SBarry Smith   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
4444f6d32bSSatish Balay /*  ISView(*rperm, STDOUT_VIEWER_SELF);*/
45f31bba2fSSatish Balay 
463d54f168SSatish Balay   /*ViewerFileOpenASCII(MPI_COMM_SELF,"row_is_orig", &V1);
47f31bba2fSSatish Balay   ViewerFileOpenASCII(MPI_COMM_SELF,"col_is_orig", &V2);
48f31bba2fSSatish Balay   ISView(*rperm,V1);
49f31bba2fSSatish Balay   ISView(*cperm,V2);
50f31bba2fSSatish Balay   ViewerDestroy(V1);
513d54f168SSatish Balay   ViewerDestroy(V2);*/
52f31bba2fSSatish Balay 
530452661fSBarry Smith   PetscFree(ia); PetscFree(ja);
5417ab2063SBarry Smith   return 0;
5517ab2063SBarry Smith }
5617ab2063SBarry Smith 
57416022c9SBarry Smith #define CHUNKSIZE   10
5817ab2063SBarry Smith 
5917ab2063SBarry Smith /* This version has row oriented v  */
60416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
6117ab2063SBarry Smith {
62416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
63416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
644b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
65d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
66416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
6717ab2063SBarry Smith 
6817ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
69416022c9SBarry Smith     row  = im[k];
7017ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
71416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
7217ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
7317ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
74416022c9SBarry Smith     low = 0;
7517ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
76416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
77416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
784b0e389bSBarry Smith       col = in[l] - shift;
794b0e389bSBarry Smith       if (roworiented) {
804b0e389bSBarry Smith         value = *v++;
814b0e389bSBarry Smith       }
824b0e389bSBarry Smith       else {
834b0e389bSBarry Smith         value = v[k + l*m];
844b0e389bSBarry Smith       }
85416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
86416022c9SBarry Smith       while (high-low > 5) {
87416022c9SBarry Smith         t = (low+high)/2;
88416022c9SBarry Smith         if (rp[t] > col) high = t;
89416022c9SBarry Smith         else             low  = t;
9017ab2063SBarry Smith       }
91416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
9217ab2063SBarry Smith         if (rp[i] > col) break;
9317ab2063SBarry Smith         if (rp[i] == col) {
94416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
9517ab2063SBarry Smith           else                  ap[i] = value;
9617ab2063SBarry Smith           goto noinsert;
9717ab2063SBarry Smith         }
9817ab2063SBarry Smith       }
9917ab2063SBarry Smith       if (nonew) goto noinsert;
10017ab2063SBarry Smith       if (nrow >= rmax) {
10117ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
102416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
10317ab2063SBarry Smith         Scalar *new_a;
10417ab2063SBarry Smith 
10517ab2063SBarry Smith         /* malloc new storage space */
106416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1070452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
10817ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
10917ab2063SBarry Smith         new_i   = new_j + new_nz;
11017ab2063SBarry Smith 
11117ab2063SBarry Smith         /* copy over old data into new slots */
11217ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
113416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
114416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
115416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
116416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
11717ab2063SBarry Smith                                                            len*sizeof(int));
118416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
119416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
12017ab2063SBarry Smith                                                            len*sizeof(Scalar));
12117ab2063SBarry Smith         /* free up old matrix storage */
1220452661fSBarry Smith         PetscFree(a->a);
1230452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
124416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
125416022c9SBarry Smith         a->singlemalloc = 1;
12617ab2063SBarry Smith 
12717ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
128416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
129416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
130416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
13117ab2063SBarry Smith       }
132416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
133416022c9SBarry Smith       /* shift up all the later entries in this row */
134416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
13517ab2063SBarry Smith         rp[ii+1] = rp[ii];
13617ab2063SBarry Smith         ap[ii+1] = ap[ii];
13717ab2063SBarry Smith       }
13817ab2063SBarry Smith       rp[i] = col;
13917ab2063SBarry Smith       ap[i] = value;
14017ab2063SBarry Smith       noinsert:;
141416022c9SBarry Smith       low = i + 1;
14217ab2063SBarry Smith     }
14317ab2063SBarry Smith     ailen[row] = nrow;
14417ab2063SBarry Smith   }
14517ab2063SBarry Smith   return 0;
14617ab2063SBarry Smith }
14717ab2063SBarry Smith 
1487eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1497eb43aa7SLois Curfman McInnes {
1507eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
151b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1527eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1537eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
1547eb43aa7SLois Curfman McInnes 
1557eb43aa7SLois Curfman McInnes   if (!a->assembled) SETERRQ(1,"MatGetValues_SeqAIJ:Not for unassembled matrix");
1567eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
1577eb43aa7SLois Curfman McInnes     row  = im[k];
1587eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
1597eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
1607eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
1617eb43aa7SLois Curfman McInnes     nrow = ailen[row];
1627eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
1637eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
1647eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
1657eb43aa7SLois Curfman McInnes       col = in[l] - shift;
1667eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
1677eb43aa7SLois Curfman McInnes       while (high-low > 5) {
1687eb43aa7SLois Curfman McInnes         t = (low+high)/2;
1697eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
1707eb43aa7SLois Curfman McInnes         else             low  = t;
1717eb43aa7SLois Curfman McInnes       }
1727eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
1737eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
1747eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
175b49de8d1SLois Curfman McInnes           *v++ = ap[i];
1767eb43aa7SLois Curfman McInnes           goto finished;
1777eb43aa7SLois Curfman McInnes         }
1787eb43aa7SLois Curfman McInnes       }
179b49de8d1SLois Curfman McInnes       *v++ = zero;
1807eb43aa7SLois Curfman McInnes       finished:;
1817eb43aa7SLois Curfman McInnes     }
1827eb43aa7SLois Curfman McInnes   }
1837eb43aa7SLois Curfman McInnes   return 0;
1847eb43aa7SLois Curfman McInnes }
1857eb43aa7SLois Curfman McInnes 
18617ab2063SBarry Smith #include "draw.h"
18717ab2063SBarry Smith #include "pinclude/pviewer.h"
188416022c9SBarry Smith #include "sysio.h"
18917ab2063SBarry Smith 
190416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
19117ab2063SBarry Smith {
192416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
193416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
19417ab2063SBarry Smith 
195416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1960452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
197416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
198416022c9SBarry Smith   col_lens[1] = a->m;
199416022c9SBarry Smith   col_lens[2] = a->n;
200416022c9SBarry Smith   col_lens[3] = a->nz;
201416022c9SBarry Smith 
202416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
203416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
204416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
20517ab2063SBarry Smith   }
206416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
2070452661fSBarry Smith   PetscFree(col_lens);
208416022c9SBarry Smith 
209416022c9SBarry Smith   /* store column indices (zero start index) */
210416022c9SBarry Smith   if (a->indexshift) {
211416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
21217ab2063SBarry Smith   }
213416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
214416022c9SBarry Smith   if (a->indexshift) {
215416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
21617ab2063SBarry Smith   }
217416022c9SBarry Smith 
218416022c9SBarry Smith   /* store nonzero values */
219416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
22017ab2063SBarry Smith   return 0;
22117ab2063SBarry Smith }
222416022c9SBarry Smith 
223416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
224416022c9SBarry Smith {
225416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
226416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
22717ab2063SBarry Smith   FILE        *fd;
22817ab2063SBarry Smith   char        *outputname;
22917ab2063SBarry Smith 
230416022c9SBarry Smith   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
231416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
232416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
23317ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
23408480c60SBarry Smith     /* no need to print additional information */ ;
23517ab2063SBarry Smith   }
23617ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
23717ab2063SBarry Smith     int nz, nzalloc, mem;
238416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
239416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
24017ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
24117ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
24217ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
24317ab2063SBarry Smith 
24417ab2063SBarry Smith     for (i=0; i<m; i++) {
245416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
24617ab2063SBarry Smith #if defined(PETSC_COMPLEX)
247416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
248416022c9SBarry Smith                    imag(a->a[j]));
24917ab2063SBarry Smith #else
250416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
25117ab2063SBarry Smith #endif
25217ab2063SBarry Smith       }
25317ab2063SBarry Smith     }
25417ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
25517ab2063SBarry Smith   }
25617ab2063SBarry Smith   else {
25717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
25817ab2063SBarry Smith       fprintf(fd,"row %d:",i);
259416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
26017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
261416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
262416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
26317ab2063SBarry Smith         }
26417ab2063SBarry Smith         else {
265416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
26617ab2063SBarry Smith         }
26717ab2063SBarry Smith #else
268416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
26917ab2063SBarry Smith #endif
27017ab2063SBarry Smith       }
27117ab2063SBarry Smith       fprintf(fd,"\n");
27217ab2063SBarry Smith     }
27317ab2063SBarry Smith   }
27417ab2063SBarry Smith   fflush(fd);
275416022c9SBarry Smith   return 0;
276416022c9SBarry Smith }
277416022c9SBarry Smith 
278416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
279416022c9SBarry Smith {
280416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
281cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
282cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
283d7e8b826SBarry Smith   Draw     draw = (Draw) viewer;
284cddf8d76SBarry Smith   DrawButton  button;
285cddf8d76SBarry Smith 
286416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
287416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
288416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
289416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
290cddf8d76SBarry Smith   color = DRAW_BLUE;
291416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
292cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
293416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
294cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
295cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
296cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
297cddf8d76SBarry Smith #else
298cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
299cddf8d76SBarry Smith #endif
300cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
301cddf8d76SBarry Smith     }
302cddf8d76SBarry Smith   }
303cddf8d76SBarry Smith   color = DRAW_CYAN;
304cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
305cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
306cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
307cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
308cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
309cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
310cddf8d76SBarry Smith     }
311cddf8d76SBarry Smith   }
312cddf8d76SBarry Smith   color = DRAW_RED;
313cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
314cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
315cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
316cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
317cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
318cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
319cddf8d76SBarry Smith #else
320cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
321cddf8d76SBarry Smith #endif
322cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
323416022c9SBarry Smith     }
324416022c9SBarry Smith   }
325416022c9SBarry Smith   DrawFlush(draw);
326cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
327cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
328cddf8d76SBarry Smith 
329cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
330cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
331cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
332cddf8d76SBarry Smith     DrawClear(draw);
333cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
334cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
335cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
336cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
337cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
338cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
339cddf8d76SBarry Smith     w *= scale; h *= scale;
340cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
341cddf8d76SBarry Smith     color = DRAW_BLUE;
342cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
343cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
344cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
345cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
346cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
347cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
348cddf8d76SBarry Smith #else
349cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
350cddf8d76SBarry Smith #endif
351cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
352cddf8d76SBarry Smith       }
353cddf8d76SBarry Smith     }
354cddf8d76SBarry Smith     color = DRAW_CYAN;
355cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
356cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
357cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
358cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
359cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
360cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
361cddf8d76SBarry Smith       }
362cddf8d76SBarry Smith     }
363cddf8d76SBarry Smith     color = DRAW_RED;
364cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
365cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
366cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
367cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
368cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
369cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
370cddf8d76SBarry Smith #else
371cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
372cddf8d76SBarry Smith #endif
373cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
374cddf8d76SBarry Smith       }
375cddf8d76SBarry Smith     }
376cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
377cddf8d76SBarry Smith   }
378416022c9SBarry Smith   return 0;
379416022c9SBarry Smith }
380416022c9SBarry Smith 
381416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
382416022c9SBarry Smith {
383416022c9SBarry Smith   Mat         A = (Mat) obj;
384416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
385416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
386416022c9SBarry Smith 
387416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
388416022c9SBarry Smith   if (!viewer) {
389416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
390416022c9SBarry Smith   }
391416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
392416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
393416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
394416022c9SBarry Smith     }
395416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
396416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
397416022c9SBarry Smith     }
398416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
399416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
400416022c9SBarry Smith     }
401416022c9SBarry Smith   }
402416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
403416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
404416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
40517ab2063SBarry Smith   }
40617ab2063SBarry Smith   return 0;
40717ab2063SBarry Smith }
40841c01911SSatish Balay int Mat_AIJ_CheckInode(Mat);
409416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
41017ab2063SBarry Smith {
411416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
41241c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
413416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
414416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
41517ab2063SBarry Smith 
41617ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
41717ab2063SBarry Smith 
41817ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
419416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
42017ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
42117ab2063SBarry Smith     if (fshift) {
422416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
42317ab2063SBarry Smith       N = ailen[i];
42417ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
42517ab2063SBarry Smith         ip[j-fshift] = ip[j];
42617ab2063SBarry Smith         ap[j-fshift] = ap[j];
42717ab2063SBarry Smith       }
42817ab2063SBarry Smith     }
42917ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
43017ab2063SBarry Smith   }
43117ab2063SBarry Smith   if (m) {
43217ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
43317ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
43417ab2063SBarry Smith   }
43517ab2063SBarry Smith   /* reset ilen and imax for each row */
43617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
43717ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
43817ab2063SBarry Smith   }
439416022c9SBarry Smith   a->nz = ai[m] + shift;
44017ab2063SBarry Smith 
44117ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
442416022c9SBarry Smith   if (fshift && a->diag) {
4430452661fSBarry Smith     PetscFree(a->diag);
444416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
445416022c9SBarry Smith     a->diag = 0;
44617ab2063SBarry Smith   }
44776dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
44841c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
449416022c9SBarry Smith   a->assembled = 1;
45017ab2063SBarry Smith   return 0;
45117ab2063SBarry Smith }
45217ab2063SBarry Smith 
453416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
45417ab2063SBarry Smith {
455416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
456cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
45717ab2063SBarry Smith   return 0;
45817ab2063SBarry Smith }
459416022c9SBarry Smith 
46017ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
46117ab2063SBarry Smith {
462416022c9SBarry Smith   Mat        A  = (Mat) obj;
463416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
464d5d45c9bSBarry Smith 
46517ab2063SBarry Smith #if defined(PETSC_LOG)
466416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
46717ab2063SBarry Smith #endif
4680452661fSBarry Smith   PetscFree(a->a);
4690452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4700452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
4710452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4720452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
4730452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
47476dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
4750452661fSBarry Smith   PetscFree(a);
476416022c9SBarry Smith   PLogObjectDestroy(A);
4770452661fSBarry Smith   PetscHeaderDestroy(A);
47817ab2063SBarry Smith   return 0;
47917ab2063SBarry Smith }
48017ab2063SBarry Smith 
481416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
48217ab2063SBarry Smith {
48317ab2063SBarry Smith   return 0;
48417ab2063SBarry Smith }
48517ab2063SBarry Smith 
486416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
48717ab2063SBarry Smith {
488416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
489416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
4904b0e389bSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
491416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
492416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
493416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
494e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
495e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
496e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
497e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
498df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
499df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
5004d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
5016c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
5026c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
5036c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
5046c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
5056c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
506e2f28af5SBarry Smith   else
5074d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
50817ab2063SBarry Smith   return 0;
50917ab2063SBarry Smith }
51017ab2063SBarry Smith 
511416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
51217ab2063SBarry Smith {
513416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
514416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
51517ab2063SBarry Smith   Scalar     *x, zero = 0.0;
51617ab2063SBarry Smith 
517416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
51817ab2063SBarry Smith   VecSet(&zero,v);
51917ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
520416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
521416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
522416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
523416022c9SBarry Smith       if (a->j[j]+shift == i) {
524416022c9SBarry Smith         x[i] = a->a[j];
52517ab2063SBarry Smith         break;
52617ab2063SBarry Smith       }
52717ab2063SBarry Smith     }
52817ab2063SBarry Smith   }
52917ab2063SBarry Smith   return 0;
53017ab2063SBarry Smith }
53117ab2063SBarry Smith 
53217ab2063SBarry Smith /* -------------------------------------------------------*/
53317ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
53417ab2063SBarry Smith /* -------------------------------------------------------*/
535416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
53617ab2063SBarry Smith {
537416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
53817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
539416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
54017ab2063SBarry Smith 
541416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
54217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
543cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
54417ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
54517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
546416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
547416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
548416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
54917ab2063SBarry Smith     alpha = x[i];
55017ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
55117ab2063SBarry Smith   }
552416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
55317ab2063SBarry Smith   return 0;
55417ab2063SBarry Smith }
555d5d45c9bSBarry Smith 
556416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
55717ab2063SBarry Smith {
558416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
55917ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
560416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
56117ab2063SBarry Smith 
562416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
56317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
56417ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
56517ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
56617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
567416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
568416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
569416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
57017ab2063SBarry Smith     alpha = x[i];
57117ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
57217ab2063SBarry Smith   }
57317ab2063SBarry Smith   return 0;
57417ab2063SBarry Smith }
57517ab2063SBarry Smith 
576416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
57717ab2063SBarry Smith {
578416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
57917ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
580416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
58117ab2063SBarry Smith 
582416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
58317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
58417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
585416022c9SBarry Smith   idx  = a->j;
586416022c9SBarry Smith   v    = a->a;
587416022c9SBarry Smith   ii   = a->i;
58817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
589416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
59017ab2063SBarry Smith     sum  = 0.0;
59117ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
59217ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
593416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
59417ab2063SBarry Smith     y[i] = sum;
59517ab2063SBarry Smith   }
596416022c9SBarry Smith   PLogFlops(2*a->nz - m);
59717ab2063SBarry Smith   return 0;
59817ab2063SBarry Smith }
59917ab2063SBarry Smith 
600416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
60117ab2063SBarry Smith {
602416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
60317ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
604cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
60517ab2063SBarry Smith 
60648d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
60717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
60817ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
609cddf8d76SBarry Smith   idx  = a->j;
610cddf8d76SBarry Smith   v    = a->a;
611cddf8d76SBarry Smith   ii   = a->i;
61217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
613cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
61417ab2063SBarry Smith     sum  = y[i];
615cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
616cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
61717ab2063SBarry Smith     z[i] = sum;
61817ab2063SBarry Smith   }
619416022c9SBarry Smith   PLogFlops(2*a->nz);
62017ab2063SBarry Smith   return 0;
62117ab2063SBarry Smith }
62217ab2063SBarry Smith 
62317ab2063SBarry Smith /*
62417ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
62517ab2063SBarry Smith */
62617ab2063SBarry Smith 
627416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
62817ab2063SBarry Smith {
629416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
630416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
63117ab2063SBarry Smith 
632416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
6330452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
634416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
635416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
636416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
637416022c9SBarry Smith       if (a->j[j]+shift == i) {
63817ab2063SBarry Smith         diag[i] = j - shift;
63917ab2063SBarry Smith         break;
64017ab2063SBarry Smith       }
64117ab2063SBarry Smith     }
64217ab2063SBarry Smith   }
643416022c9SBarry Smith   a->diag = diag;
64417ab2063SBarry Smith   return 0;
64517ab2063SBarry Smith }
64617ab2063SBarry Smith 
647416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
64817ab2063SBarry Smith                            double fshift,int its,Vec xx)
64917ab2063SBarry Smith {
650416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
651416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
652d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
65317ab2063SBarry Smith 
65417ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
655416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
656416022c9SBarry Smith   diag = a->diag;
657416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
65817ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
65917ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
66017ab2063SBarry Smith     bs = b + shift;
66117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
662416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
663416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
664416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
665416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
66617ab2063SBarry Smith         sum  = b[i]*d/omega;
66717ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
66817ab2063SBarry Smith         x[i] = sum;
66917ab2063SBarry Smith     }
67017ab2063SBarry Smith     return 0;
67117ab2063SBarry Smith   }
67217ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
67317ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
67417ab2063SBarry Smith   }
675416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
67617ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
67717ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
67817ab2063SBarry Smith 
67917ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
68017ab2063SBarry Smith 
68117ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
68217ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
68317ab2063SBarry Smith     is the relaxation factor.
68417ab2063SBarry Smith     */
6850452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
68617ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
68717ab2063SBarry Smith 
68817ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
68917ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
690416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
691416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
692416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
693416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
69417ab2063SBarry Smith       sum  = b[i];
69517ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
69617ab2063SBarry Smith       x[i] = omega*(sum/d);
69717ab2063SBarry Smith     }
69817ab2063SBarry Smith 
69917ab2063SBarry Smith     /*  t = b - (2*E - D)x */
700416022c9SBarry Smith     v = a->a;
70117ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
70217ab2063SBarry Smith 
70317ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
704416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
705416022c9SBarry Smith     diag = a->diag;
70617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
707416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
708416022c9SBarry Smith       n    = diag[i] - a->i[i];
709416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
710416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
71117ab2063SBarry Smith       sum  = t[i];
71217ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
71317ab2063SBarry Smith       t[i] = omega*(sum/d);
71417ab2063SBarry Smith     }
71517ab2063SBarry Smith 
71617ab2063SBarry Smith     /*  x = x + t */
71717ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7180452661fSBarry Smith     PetscFree(t);
71917ab2063SBarry Smith     return 0;
72017ab2063SBarry Smith   }
72117ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
72217ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
72317ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
724416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
725416022c9SBarry Smith         n    = diag[i] - a->i[i];
726416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
727416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
72817ab2063SBarry Smith         sum  = b[i];
72917ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
73017ab2063SBarry Smith         x[i] = omega*(sum/d);
73117ab2063SBarry Smith       }
73217ab2063SBarry Smith       xb = x;
73317ab2063SBarry Smith     }
73417ab2063SBarry Smith     else xb = b;
73517ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
73617ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
73717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
738416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
73917ab2063SBarry Smith       }
74017ab2063SBarry Smith     }
74117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
74217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
743416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
744416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
745416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
746416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
74717ab2063SBarry Smith         sum  = xb[i];
74817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
74917ab2063SBarry Smith         x[i] = omega*(sum/d);
75017ab2063SBarry Smith       }
75117ab2063SBarry Smith     }
75217ab2063SBarry Smith     its--;
75317ab2063SBarry Smith   }
75417ab2063SBarry Smith   while (its--) {
75517ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
75617ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
757416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
758416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
759416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
760416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
76117ab2063SBarry Smith         sum  = b[i];
76217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
76317ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
76417ab2063SBarry Smith       }
76517ab2063SBarry Smith     }
76617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
76717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
768416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
769416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
770416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
771416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
77217ab2063SBarry Smith         sum  = b[i];
77317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
77417ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
77517ab2063SBarry Smith       }
77617ab2063SBarry Smith     }
77717ab2063SBarry Smith   }
77817ab2063SBarry Smith   return 0;
77917ab2063SBarry Smith }
78017ab2063SBarry Smith 
781d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
78217ab2063SBarry Smith {
783416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
784416022c9SBarry Smith   *nz      = a->nz;
785416022c9SBarry Smith   *nzalloc = a->maxnz;
786416022c9SBarry Smith   *mem     = (int)A->mem;
78717ab2063SBarry Smith   return 0;
78817ab2063SBarry Smith }
78917ab2063SBarry Smith 
79017ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
79117ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
79217ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
79317ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
79417ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
79517ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
79617ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
79717ab2063SBarry Smith 
79817ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
79917ab2063SBarry Smith {
800416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
801416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
80217ab2063SBarry Smith 
80317ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
80417ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
80517ab2063SBarry Smith   if (diag) {
80617ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
807416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
808416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
809416022c9SBarry Smith         a->ilen[rows[i]] = 1;
810416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
811416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
81217ab2063SBarry Smith       }
81317ab2063SBarry Smith       else {
81417ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
81517ab2063SBarry Smith         CHKERRQ(ierr);
81617ab2063SBarry Smith       }
81717ab2063SBarry Smith     }
81817ab2063SBarry Smith   }
81917ab2063SBarry Smith   else {
82017ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
821416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
822416022c9SBarry Smith       a->ilen[rows[i]] = 0;
82317ab2063SBarry Smith     }
82417ab2063SBarry Smith   }
82517ab2063SBarry Smith   ISRestoreIndices(is,&rows);
82617ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82717ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82817ab2063SBarry Smith   return 0;
82917ab2063SBarry Smith }
83017ab2063SBarry Smith 
831416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
83217ab2063SBarry Smith {
833416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
834416022c9SBarry Smith   *m = a->m; *n = a->n;
83517ab2063SBarry Smith   return 0;
83617ab2063SBarry Smith }
83717ab2063SBarry Smith 
838416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
83917ab2063SBarry Smith {
840416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
841416022c9SBarry Smith   *m = 0; *n = a->m;
84217ab2063SBarry Smith   return 0;
84317ab2063SBarry Smith }
844416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
84517ab2063SBarry Smith {
846416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
847416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
84817ab2063SBarry Smith 
849416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
85017ab2063SBarry Smith 
851416022c9SBarry Smith   if (!a->assembled) {
852416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
853416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
85417ab2063SBarry Smith   }
855416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
856416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
85717ab2063SBarry Smith   if (idx) {
85817ab2063SBarry Smith     if (*nz) {
859416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
8600452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
86117ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
86217ab2063SBarry Smith     }
86317ab2063SBarry Smith     else *idx = 0;
86417ab2063SBarry Smith   }
86517ab2063SBarry Smith   return 0;
86617ab2063SBarry Smith }
86717ab2063SBarry Smith 
868416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
86917ab2063SBarry Smith {
8700452661fSBarry Smith   if (idx) {if (*idx) PetscFree(*idx);}
87117ab2063SBarry Smith   return 0;
87217ab2063SBarry Smith }
87317ab2063SBarry Smith 
874cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
87517ab2063SBarry Smith {
876416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
877416022c9SBarry Smith   Scalar     *v = a->a;
87817ab2063SBarry Smith   double     sum = 0.0;
879416022c9SBarry Smith   int        i, j,shift = a->indexshift;
88017ab2063SBarry Smith 
881416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
88217ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
883416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
88417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
88517ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
88617ab2063SBarry Smith #else
88717ab2063SBarry Smith       sum += (*v)*(*v); v++;
88817ab2063SBarry Smith #endif
88917ab2063SBarry Smith     }
89017ab2063SBarry Smith     *norm = sqrt(sum);
89117ab2063SBarry Smith   }
89217ab2063SBarry Smith   else if (type == NORM_1) {
89317ab2063SBarry Smith     double *tmp;
894416022c9SBarry Smith     int    *jj = a->j;
8950452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
896cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
89717ab2063SBarry Smith     *norm = 0.0;
898416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
899a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
90017ab2063SBarry Smith     }
901416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
90217ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
90317ab2063SBarry Smith     }
9040452661fSBarry Smith     PetscFree(tmp);
90517ab2063SBarry Smith   }
90617ab2063SBarry Smith   else if (type == NORM_INFINITY) {
90717ab2063SBarry Smith     *norm = 0.0;
908416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
909416022c9SBarry Smith       v = a->a + a->i[j] + shift;
91017ab2063SBarry Smith       sum = 0.0;
911416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
912cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
91317ab2063SBarry Smith       }
91417ab2063SBarry Smith       if (sum > *norm) *norm = sum;
91517ab2063SBarry Smith     }
91617ab2063SBarry Smith   }
91717ab2063SBarry Smith   else {
91848d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
91917ab2063SBarry Smith   }
92017ab2063SBarry Smith   return 0;
92117ab2063SBarry Smith }
92217ab2063SBarry Smith 
923416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
92417ab2063SBarry Smith {
925416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
926416022c9SBarry Smith   Mat        C;
927416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
928416022c9SBarry Smith   int        shift = a->indexshift;
929d5d45c9bSBarry Smith   Scalar     *array = a->a;
93017ab2063SBarry Smith 
931416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
9320452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
933cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
93417ab2063SBarry Smith   if (shift) {
93517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
93617ab2063SBarry Smith   }
93717ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
938416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
9390452661fSBarry Smith   PetscFree(col);
94017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
94117ab2063SBarry Smith     len = ai[i+1]-ai[i];
942416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
94317ab2063SBarry Smith     array += len; aj += len;
94417ab2063SBarry Smith   }
94517ab2063SBarry Smith   if (shift) {
94617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
94717ab2063SBarry Smith   }
94817ab2063SBarry Smith 
949416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
950416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
95117ab2063SBarry Smith 
952416022c9SBarry Smith   if (B) {
953416022c9SBarry Smith     *B = C;
95417ab2063SBarry Smith   } else {
955416022c9SBarry Smith     /* This isn't really an in-place transpose */
9560452661fSBarry Smith     PetscFree(a->a);
9570452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
9580452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
9590452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
9600452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
9610452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
9620452661fSBarry Smith     PetscFree(a);
963416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
9640452661fSBarry Smith     PetscHeaderDestroy(C);
96517ab2063SBarry Smith   }
96617ab2063SBarry Smith   return 0;
96717ab2063SBarry Smith }
96817ab2063SBarry Smith 
969416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
97017ab2063SBarry Smith {
971416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
97217ab2063SBarry Smith   Scalar     *l,*r,x,*v;
973d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
97417ab2063SBarry Smith 
97548d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
97617ab2063SBarry Smith   if (ll) {
97717ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
978416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
979416022c9SBarry Smith     v = a->a;
98017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
98117ab2063SBarry Smith       x = l[i];
982416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
98317ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
98417ab2063SBarry Smith     }
98517ab2063SBarry Smith   }
98617ab2063SBarry Smith   if (rr) {
98717ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
988416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
989416022c9SBarry Smith     v = a->a; jj = a->j;
99017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
99117ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
99217ab2063SBarry Smith     }
99317ab2063SBarry Smith   }
99417ab2063SBarry Smith   return 0;
99517ab2063SBarry Smith }
99617ab2063SBarry Smith 
997cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
99817ab2063SBarry Smith {
999db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
100002834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
1001a2744918SBarry Smith   register int sum,lensi;
100202834360SBarry Smith   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
1003a2744918SBarry Smith   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
1004db02288aSLois Curfman McInnes   Scalar       *vwork,*a_new;
1005416022c9SBarry Smith   Mat          C;
100617ab2063SBarry Smith 
1007416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
100817ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
100917ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
101017ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
101117ab2063SBarry Smith 
101202834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
101302834360SBarry Smith     /* special case of contiguous rows */
10140452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
101502834360SBarry Smith     starts = lens + ncols;
101602834360SBarry Smith     /* loop over new rows determining lens and starting points */
101702834360SBarry Smith     for (i=0; i<nrows; i++) {
1018a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1019a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
102002834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1021*d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
102202834360SBarry Smith           starts[i] = k;
102302834360SBarry Smith           break;
102402834360SBarry Smith 	}
102502834360SBarry Smith       }
1026a2744918SBarry Smith       sum = 0;
102702834360SBarry Smith       while (k < kend) {
1028*d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1029a2744918SBarry Smith         sum++;
103002834360SBarry Smith       }
1031a2744918SBarry Smith       lens[i] = sum;
103202834360SBarry Smith     }
103302834360SBarry Smith     /* create submatrix */
1034cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
103508480c60SBarry Smith       int n_cols,n_rows;
103608480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
103708480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1038*d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
103908480c60SBarry Smith       C = *B;
104008480c60SBarry Smith     }
104108480c60SBarry Smith     else {
104202834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
104308480c60SBarry Smith     }
1044db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1045db02288aSLois Curfman McInnes 
104602834360SBarry Smith     /* loop over rows inserting into submatrix */
1047db02288aSLois Curfman McInnes     a_new    = c->a;
1048db02288aSLois Curfman McInnes     j_new    = c->j;
1049db02288aSLois Curfman McInnes     i_new    = c->i;
1050db02288aSLois Curfman McInnes     i_new[0] = -shift;
105102834360SBarry Smith     for (i=0; i<nrows; i++) {
1052a2744918SBarry Smith       ii    = starts[i];
1053a2744918SBarry Smith       lensi = lens[i];
1054a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1055a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
105602834360SBarry Smith       }
1057a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1058a2744918SBarry Smith       a_new      += lensi;
1059a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1060a2744918SBarry Smith       c->ilen[i]  = lensi;
106102834360SBarry Smith     }
10620452661fSBarry Smith     PetscFree(lens);
106302834360SBarry Smith   }
106402834360SBarry Smith   else {
106502834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10660452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
106702834360SBarry Smith     ssmap = smap + shift;
10687eb43aa7SLois Curfman McInnes     cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork);
106902834360SBarry Smith     lens  = cwork + ncols;
10700452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1071cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
107217ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
107302834360SBarry Smith     /* determine lens of each row */
107402834360SBarry Smith     for (i=0; i<nrows; i++) {
1075*d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
107602834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
107702834360SBarry Smith       lens[i] = 0;
107802834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1079*d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
108002834360SBarry Smith           lens[i]++;
108102834360SBarry Smith         }
108202834360SBarry Smith       }
108302834360SBarry Smith     }
108417ab2063SBarry Smith     /* Create and fill new matrix */
1085a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
108608480c60SBarry Smith       int n_cols,n_rows;
108708480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
108808480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1089*d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
109008480c60SBarry Smith       C = *B;
109108480c60SBarry Smith     }
109208480c60SBarry Smith     else {
109302834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
109408480c60SBarry Smith     }
109517ab2063SBarry Smith     for (i=0; i<nrows; i++) {
109617ab2063SBarry Smith       nznew  = 0;
1097*d8ced48eSBarry Smith       kstart = ai[irow[i]]+shift;
1098416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
109917ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
110002834360SBarry Smith         if (ssmap[a->j[k]]) {
110102834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1102416022c9SBarry Smith           vwork[nznew++] = a->a[k];
110317ab2063SBarry Smith         }
110417ab2063SBarry Smith       }
1105416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
110617ab2063SBarry Smith     }
110702834360SBarry Smith     /* Free work space */
110802834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
11090452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
111002834360SBarry Smith   }
1111416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1112416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
111317ab2063SBarry Smith 
111417ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1115416022c9SBarry Smith   *B = C;
111617ab2063SBarry Smith   return 0;
111717ab2063SBarry Smith }
111817ab2063SBarry Smith 
1119a871dcd8SBarry Smith /*
112063b91edcSBarry Smith      note: This can only work for identity for row and col. It would
112163b91edcSBarry Smith    be good to check this and otherwise generate an error.
1122a871dcd8SBarry Smith */
112363b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1124a871dcd8SBarry Smith {
112563b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
112608480c60SBarry Smith   int        ierr;
112763b91edcSBarry Smith   Mat        outA;
112863b91edcSBarry Smith 
1129a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1130a871dcd8SBarry Smith 
113163b91edcSBarry Smith   outA          = inA;
113263b91edcSBarry Smith   inA->factor   = FACTOR_LU;
113363b91edcSBarry Smith   a->row        = row;
113463b91edcSBarry Smith   a->col        = col;
113563b91edcSBarry Smith 
11360452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
113763b91edcSBarry Smith 
113808480c60SBarry Smith   if (!a->diag) {
113908480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
114063b91edcSBarry Smith   }
114163b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1142a871dcd8SBarry Smith   return 0;
1143a871dcd8SBarry Smith }
1144a871dcd8SBarry Smith 
1145cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1146cddf8d76SBarry Smith                                     Mat **B)
1147cddf8d76SBarry Smith {
1148cddf8d76SBarry Smith   int ierr,i;
1149cddf8d76SBarry Smith 
1150cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
11510452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1152cddf8d76SBarry Smith   }
1153cddf8d76SBarry Smith 
1154cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1155cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1156cddf8d76SBarry Smith   }
1157cddf8d76SBarry Smith   return 0;
1158cddf8d76SBarry Smith }
1159cddf8d76SBarry Smith 
11604dcbc457SBarry Smith static int MatIncreaseOverlap_SeqAIJ(Mat A, int n, IS *is, int ov)
11614dcbc457SBarry Smith {
11624dcbc457SBarry Smith   int i,m,*idx,ierr;
11634dcbc457SBarry Smith 
11644dcbc457SBarry Smith   for ( i=0; i<n; i++ ) {
11654dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr);
11664dcbc457SBarry Smith     ISGetLocalSize(is[i],&m);
11674dcbc457SBarry Smith   }
11684dcbc457SBarry Smith   SETERRQ(1,"MatIncreaseOverlap_SeqAIJ:Not implemented");
11694dcbc457SBarry Smith }
117017ab2063SBarry Smith 
1171682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1172682d7d0cSBarry Smith {
1173682d7d0cSBarry Smith   static int called = 0;
1174682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1175682d7d0cSBarry Smith 
1176682d7d0cSBarry Smith   if (called) return 0; else called = 1;
11772a7368beSLois Curfman McInnes   MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1178682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
11792a7368beSLois Curfman McInnes   MPIU_printf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
1180682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
1181682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1182682d7d0cSBarry Smith #if defined(HAVE_ESSL)
1183682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1184682d7d0cSBarry Smith #endif
1185682d7d0cSBarry Smith   return 0;
1186682d7d0cSBarry Smith }
1187682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
118817ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
118917ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1190416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1191416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
119217ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
119317ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
119417ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
119517ab2063SBarry Smith        MatRelax_SeqAIJ,
119617ab2063SBarry Smith        MatTranspose_SeqAIJ,
119717ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
119817ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
119917ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
120017ab2063SBarry Smith        MatCompress_SeqAIJ,
120117ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
120217ab2063SBarry Smith        MatGetReordering_SeqAIJ,
120317ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
120417ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
120517ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
120617ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1207416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
12083d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1209cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
12107eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1211682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1212682d7d0cSBarry Smith        MatPrintHelp_SeqAIJ};
121317ab2063SBarry Smith 
121417ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
121517ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
121617ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
121717ab2063SBarry Smith 
121817ab2063SBarry Smith /*@C
1219682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
122017ab2063SBarry Smith                      (the default uniprocessor PETSc format).
122117ab2063SBarry Smith 
122217ab2063SBarry Smith    Input Parameters:
122317ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
122417ab2063SBarry Smith .  m - number of rows
122517ab2063SBarry Smith .  n - number of columns
122617ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
122717ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
122817ab2063SBarry Smith 
122917ab2063SBarry Smith    Output Parameter:
1230416022c9SBarry Smith .  A - the matrix
123117ab2063SBarry Smith 
123217ab2063SBarry Smith    Notes:
123317ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
123417ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
12350002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
12360002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
123717ab2063SBarry Smith 
123817ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1239a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
124017ab2063SBarry Smith    allocation.
124117ab2063SBarry Smith 
1242682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1243682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1244682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
12456c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
12466c7ebb05SLois Curfman McInnes 
12476c7ebb05SLois Curfman McInnes    Options Database Keys:
12486c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
12496c7ebb05SLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)
12504b0e389bSBarry Smith $    -mat_aij_oneindex - Internally using indexing starting at 1 rather then zero.
12514b0e389bSBarry Smith .                        Note: You still index entries starting at 0!
125217ab2063SBarry Smith 
125317ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
125417ab2063SBarry Smith @*/
1255416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
125617ab2063SBarry Smith {
1257416022c9SBarry Smith   Mat        B;
1258416022c9SBarry Smith   Mat_SeqAIJ *b;
125917ab2063SBarry Smith   int        i,len,ierr;
1260d5d45c9bSBarry Smith 
1261416022c9SBarry Smith   *A      = 0;
12620452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1263416022c9SBarry Smith   PLogObjectCreate(B);
12640452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1265416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1266416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1267416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1268416022c9SBarry Smith   B->factor           = 0;
1269416022c9SBarry Smith   B->lupivotthreshold = 1.0;
1270b4fd4287SBarry Smith   OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1271416022c9SBarry Smith   b->row              = 0;
1272416022c9SBarry Smith   b->col              = 0;
1273416022c9SBarry Smith   b->indexshift       = 0;
1274b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_oneindex")) b->indexshift = -1;
127517ab2063SBarry Smith 
1276416022c9SBarry Smith   b->m       = m;
1277416022c9SBarry Smith   b->n       = n;
12780452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1279b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
12807b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
12817b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1282416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
128317ab2063SBarry Smith     nz = nz*m;
128417ab2063SBarry Smith   }
128517ab2063SBarry Smith   else {
128617ab2063SBarry Smith     nz = 0;
1287416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
128817ab2063SBarry Smith   }
128917ab2063SBarry Smith 
129017ab2063SBarry Smith   /* allocate the matrix space */
1291416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
12920452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1293416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1294cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1295416022c9SBarry Smith   b->i  = b->j + nz;
1296416022c9SBarry Smith   b->singlemalloc = 1;
129717ab2063SBarry Smith 
1298416022c9SBarry Smith   b->i[0] = -b->indexshift;
129917ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1300416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
130117ab2063SBarry Smith   }
130217ab2063SBarry Smith 
1303416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
13040452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1305416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1306416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
130717ab2063SBarry Smith 
1308416022c9SBarry Smith   b->nz               = 0;
1309416022c9SBarry Smith   b->maxnz            = nz;
1310416022c9SBarry Smith   b->sorted           = 0;
1311416022c9SBarry Smith   b->roworiented      = 1;
1312416022c9SBarry Smith   b->nonew            = 0;
1313416022c9SBarry Smith   b->diag             = 0;
1314416022c9SBarry Smith   b->assembled        = 0;
1315416022c9SBarry Smith   b->solve_work       = 0;
131671bd300dSLois Curfman McInnes   b->spptr            = 0;
1317754ec7b1SSatish Balay   b->inode.node_count = 0;
1318754ec7b1SSatish Balay   b->inode.size       = 0;
13196c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
13206c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
132117ab2063SBarry Smith 
1322416022c9SBarry Smith   *A = B;
1323b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_superlu")) {
1324416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
132517ab2063SBarry Smith   }
1326b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_essl")) {
1327416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
132817ab2063SBarry Smith   }
1329b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_dxml")) {
1330416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1331416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
133217ab2063SBarry Smith   }
1333682d7d0cSBarry Smith   if (OptionsHasName(PETSC_NULL,"-help")) {
1334682d7d0cSBarry Smith     ierr = MatPrintHelp(B); CHKERRQ(ierr);
1335682d7d0cSBarry Smith   }
133617ab2063SBarry Smith   return 0;
133717ab2063SBarry Smith }
133817ab2063SBarry Smith 
13393d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
134017ab2063SBarry Smith {
1341416022c9SBarry Smith   Mat        C;
1342416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
134308480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
134417ab2063SBarry Smith 
13454043dd9cSLois Curfman McInnes   *B = 0;
13463d1612f7SBarry Smith   if (!a->assembled) SETERRQ(1,"MatConvertSameType_SeqAIJ:Cannot copy unassembled matrix");
13470452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1348416022c9SBarry Smith   PLogObjectCreate(C);
13490452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
135041c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1351416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1352416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1353416022c9SBarry Smith   C->factor     = A->factor;
1354416022c9SBarry Smith   c->row        = 0;
1355416022c9SBarry Smith   c->col        = 0;
1356416022c9SBarry Smith   c->indexshift = shift;
135717ab2063SBarry Smith 
1358416022c9SBarry Smith   c->m          = a->m;
1359416022c9SBarry Smith   c->n          = a->n;
136017ab2063SBarry Smith 
13610452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
13620452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
136317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1364416022c9SBarry Smith     c->imax[i] = a->imax[i];
1365416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
136617ab2063SBarry Smith   }
136717ab2063SBarry Smith 
136817ab2063SBarry Smith   /* allocate the matrix space */
1369416022c9SBarry Smith   c->singlemalloc = 1;
1370416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
13710452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1372416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1373416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1374416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
137517ab2063SBarry Smith   if (m > 0) {
1376416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
137708480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1378416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
137917ab2063SBarry Smith     }
138008480c60SBarry Smith   }
138117ab2063SBarry Smith 
1382416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1383416022c9SBarry Smith   c->sorted      = a->sorted;
1384416022c9SBarry Smith   c->roworiented = a->roworiented;
1385416022c9SBarry Smith   c->nonew       = a->nonew;
138617ab2063SBarry Smith 
1387416022c9SBarry Smith   if (a->diag) {
13880452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1389416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
139017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1391416022c9SBarry Smith       c->diag[i] = a->diag[i];
139217ab2063SBarry Smith     }
139317ab2063SBarry Smith   }
1394416022c9SBarry Smith   else c->diag          = 0;
13956c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
13966c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1397754ec7b1SSatish Balay   if (a->inode.size){
1398754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1399754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1400754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1401754ec7b1SSatish Balay   } else {
1402754ec7b1SSatish Balay     c->inode.size       = 0;
1403754ec7b1SSatish Balay     c->inode.node_count = 0;
1404754ec7b1SSatish Balay   }
1405416022c9SBarry Smith   c->assembled          = 1;
1406416022c9SBarry Smith   c->nz                 = a->nz;
1407416022c9SBarry Smith   c->maxnz              = a->maxnz;
1408416022c9SBarry Smith   c->solve_work         = 0;
140976dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1410754ec7b1SSatish Balay 
1411416022c9SBarry Smith   *B = C;
141217ab2063SBarry Smith   return 0;
141317ab2063SBarry Smith }
141417ab2063SBarry Smith 
1415416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
141617ab2063SBarry Smith {
1417416022c9SBarry Smith   Mat_SeqAIJ   *a;
1418416022c9SBarry Smith   Mat          B;
141917699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
142017ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
142117ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
142217ab2063SBarry Smith 
142317699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
142417699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
142517ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1426416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
142748d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
142817ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
142917ab2063SBarry Smith 
143017ab2063SBarry Smith   /* read in row lengths */
14310452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1432416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
143317ab2063SBarry Smith 
143417ab2063SBarry Smith   /* create our matrix */
1435416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1436416022c9SBarry Smith   B = *A;
1437416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1438416022c9SBarry Smith   shift = a->indexshift;
143917ab2063SBarry Smith 
144017ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1441416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
144217ab2063SBarry Smith   if (shift) {
144317ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1444416022c9SBarry Smith       a->j[i] += 1;
144517ab2063SBarry Smith     }
144617ab2063SBarry Smith   }
144717ab2063SBarry Smith 
144817ab2063SBarry Smith   /* read in nonzero values */
1449416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
145017ab2063SBarry Smith 
145117ab2063SBarry Smith   /* set matrix "i" values */
1452416022c9SBarry Smith   a->i[0] = -shift;
145317ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1454416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1455416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
145617ab2063SBarry Smith   }
14570452661fSBarry Smith   PetscFree(rowlengths);
145817ab2063SBarry Smith 
1459416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1460416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
146117ab2063SBarry Smith   return 0;
146217ab2063SBarry Smith }
146317ab2063SBarry Smith 
146417ab2063SBarry Smith 
146517ab2063SBarry Smith 
1466