xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 4b14c69e1e8458494734dcd7bdc4eefb1c50097c)
117ab2063SBarry Smith #ifndef lint
2*4b14c69eSBarry Smith static char vcid[] = "$Id: aij.c,v 1.140 1996/01/26 04:33:46 bsmith 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"
12e4d965acSSatish Balay #include "petsc.h"
1317ab2063SBarry Smith 
1417ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**);
1517ab2063SBarry Smith 
16416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
1717ab2063SBarry Smith {
18416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
19a2744918SBarry Smith   int        ierr, *ia, *ja,n,*idx,i;
203d54f168SSatish Balay   /*Viewer     V1, V2;*/
2117ab2063SBarry Smith 
22a2744918SBarry Smith   /*
23a2744918SBarry Smith      this is tacky: In the future when we have written special factorization
24a2744918SBarry Smith      and solve routines for the identity permutation we should use a
25a2744918SBarry Smith      stride index set instead of the general one.
26a2744918SBarry Smith   */
27a2744918SBarry Smith   if (type  == ORDER_NATURAL) {
28a2744918SBarry Smith     n = a->n;
29a2744918SBarry Smith     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
30a2744918SBarry Smith     for ( i=0; i<n; i++ ) idx[i] = i;
31a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
32a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
33a2744918SBarry Smith     PetscFree(idx);
34a2744918SBarry Smith     ISSetPermutation(*rperm);
35a2744918SBarry Smith     ISSetPermutation(*cperm);
36a2744918SBarry Smith     ISSetIdentity(*rperm);
37a2744918SBarry Smith     ISSetIdentity(*cperm);
38a2744918SBarry Smith     return 0;
39a2744918SBarry Smith   }
40a2744918SBarry Smith 
41416022c9SBarry Smith   ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr);
42416022c9SBarry Smith   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
43f31bba2fSSatish Balay 
440452661fSBarry Smith   PetscFree(ia); PetscFree(ja);
4517ab2063SBarry Smith   return 0;
4617ab2063SBarry Smith }
4717ab2063SBarry Smith 
48227d817aSBarry Smith #define CHUNKSIZE   15
4917ab2063SBarry Smith 
5017ab2063SBarry Smith /* This version has row oriented v  */
51416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
5217ab2063SBarry Smith {
53416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
54416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
554b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
56d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
57416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
5817ab2063SBarry Smith 
5917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
60416022c9SBarry Smith     row  = im[k];
6117ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
62416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
6317ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
6417ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
65416022c9SBarry Smith     low = 0;
6617ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
67416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
68416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
694b0e389bSBarry Smith       col = in[l] - shift;
704b0e389bSBarry Smith       if (roworiented) {
714b0e389bSBarry Smith         value = *v++;
724b0e389bSBarry Smith       }
734b0e389bSBarry Smith       else {
744b0e389bSBarry Smith         value = v[k + l*m];
754b0e389bSBarry Smith       }
76416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
77416022c9SBarry Smith       while (high-low > 5) {
78416022c9SBarry Smith         t = (low+high)/2;
79416022c9SBarry Smith         if (rp[t] > col) high = t;
80416022c9SBarry Smith         else             low  = t;
8117ab2063SBarry Smith       }
82416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
8317ab2063SBarry Smith         if (rp[i] > col) break;
8417ab2063SBarry Smith         if (rp[i] == col) {
85416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
8617ab2063SBarry Smith           else                  ap[i] = value;
8717ab2063SBarry Smith           goto noinsert;
8817ab2063SBarry Smith         }
8917ab2063SBarry Smith       }
9017ab2063SBarry Smith       if (nonew) goto noinsert;
9117ab2063SBarry Smith       if (nrow >= rmax) {
9217ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
93416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
9417ab2063SBarry Smith         Scalar *new_a;
9517ab2063SBarry Smith 
9617ab2063SBarry Smith         /* malloc new storage space */
97416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
980452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
9917ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
10017ab2063SBarry Smith         new_i   = new_j + new_nz;
10117ab2063SBarry Smith 
10217ab2063SBarry Smith         /* copy over old data into new slots */
10317ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
104416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
105416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
106416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
107416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
10817ab2063SBarry Smith                                                            len*sizeof(int));
109416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
110416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
11117ab2063SBarry Smith                                                            len*sizeof(Scalar));
11217ab2063SBarry Smith         /* free up old matrix storage */
1130452661fSBarry Smith         PetscFree(a->a);
1140452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
115416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
116416022c9SBarry Smith         a->singlemalloc = 1;
11717ab2063SBarry Smith 
11817ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
119416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
120416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
121416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
12217ab2063SBarry Smith       }
123416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
124416022c9SBarry Smith       /* shift up all the later entries in this row */
125416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
12617ab2063SBarry Smith         rp[ii+1] = rp[ii];
12717ab2063SBarry Smith         ap[ii+1] = ap[ii];
12817ab2063SBarry Smith       }
12917ab2063SBarry Smith       rp[i] = col;
13017ab2063SBarry Smith       ap[i] = value;
13117ab2063SBarry Smith       noinsert:;
132416022c9SBarry Smith       low = i + 1;
13317ab2063SBarry Smith     }
13417ab2063SBarry Smith     ailen[row] = nrow;
13517ab2063SBarry Smith   }
13617ab2063SBarry Smith   return 0;
13717ab2063SBarry Smith }
13817ab2063SBarry Smith 
1397eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1407eb43aa7SLois Curfman McInnes {
1417eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
142b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1437eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1447eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
1457eb43aa7SLois Curfman McInnes 
1467eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
1477eb43aa7SLois Curfman McInnes     row  = im[k];
1487eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
1497eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
1507eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
1517eb43aa7SLois Curfman McInnes     nrow = ailen[row];
1527eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
1537eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
1547eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
1557eb43aa7SLois Curfman McInnes       col = in[l] - shift;
1567eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
1577eb43aa7SLois Curfman McInnes       while (high-low > 5) {
1587eb43aa7SLois Curfman McInnes         t = (low+high)/2;
1597eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
1607eb43aa7SLois Curfman McInnes         else             low  = t;
1617eb43aa7SLois Curfman McInnes       }
1627eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
1637eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
1647eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
165b49de8d1SLois Curfman McInnes           *v++ = ap[i];
1667eb43aa7SLois Curfman McInnes           goto finished;
1677eb43aa7SLois Curfman McInnes         }
1687eb43aa7SLois Curfman McInnes       }
169b49de8d1SLois Curfman McInnes       *v++ = zero;
1707eb43aa7SLois Curfman McInnes       finished:;
1717eb43aa7SLois Curfman McInnes     }
1727eb43aa7SLois Curfman McInnes   }
1737eb43aa7SLois Curfman McInnes   return 0;
1747eb43aa7SLois Curfman McInnes }
1757eb43aa7SLois Curfman McInnes 
17617ab2063SBarry Smith #include "draw.h"
17717ab2063SBarry Smith #include "pinclude/pviewer.h"
178416022c9SBarry Smith #include "sysio.h"
17917ab2063SBarry Smith 
180416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
18117ab2063SBarry Smith {
182416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
183416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
18417ab2063SBarry Smith 
185416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1860452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
187416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
188416022c9SBarry Smith   col_lens[1] = a->m;
189416022c9SBarry Smith   col_lens[2] = a->n;
190416022c9SBarry Smith   col_lens[3] = a->nz;
191416022c9SBarry Smith 
192416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
193416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
194416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
19517ab2063SBarry Smith   }
196416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
1970452661fSBarry Smith   PetscFree(col_lens);
198416022c9SBarry Smith 
199416022c9SBarry Smith   /* store column indices (zero start index) */
200416022c9SBarry Smith   if (a->indexshift) {
201416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
20217ab2063SBarry Smith   }
203416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
204416022c9SBarry Smith   if (a->indexshift) {
205416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
20617ab2063SBarry Smith   }
207416022c9SBarry Smith 
208416022c9SBarry Smith   /* store nonzero values */
209416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
21017ab2063SBarry Smith   return 0;
21117ab2063SBarry Smith }
212416022c9SBarry Smith 
213416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
214416022c9SBarry Smith {
215416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
216416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
21717ab2063SBarry Smith   FILE        *fd;
21817ab2063SBarry Smith   char        *outputname;
21917ab2063SBarry Smith 
220227d817aSBarry Smith   ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr);
221416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
222416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
22317ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
22408480c60SBarry Smith     /* no need to print additional information */ ;
22517ab2063SBarry Smith   }
22617ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
22717ab2063SBarry Smith     int nz, nzalloc, mem;
228416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
229416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
23017ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
23117ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
23217ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
23317ab2063SBarry Smith 
23417ab2063SBarry Smith     for (i=0; i<m; i++) {
235416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
23617ab2063SBarry Smith #if defined(PETSC_COMPLEX)
237416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
238416022c9SBarry Smith                    imag(a->a[j]));
23917ab2063SBarry Smith #else
240416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
24117ab2063SBarry Smith #endif
24217ab2063SBarry Smith       }
24317ab2063SBarry Smith     }
24417ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
24517ab2063SBarry Smith   }
24617ab2063SBarry Smith   else {
24717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
24817ab2063SBarry Smith       fprintf(fd,"row %d:",i);
249416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
25017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
251416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
252416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
25317ab2063SBarry Smith         }
25417ab2063SBarry Smith         else {
255416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
25617ab2063SBarry Smith         }
25717ab2063SBarry Smith #else
258416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
25917ab2063SBarry Smith #endif
26017ab2063SBarry Smith       }
26117ab2063SBarry Smith       fprintf(fd,"\n");
26217ab2063SBarry Smith     }
26317ab2063SBarry Smith   }
26417ab2063SBarry Smith   fflush(fd);
265416022c9SBarry Smith   return 0;
266416022c9SBarry Smith }
267416022c9SBarry Smith 
268416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
269416022c9SBarry Smith {
270416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
271cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
272cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
273d7e8b826SBarry Smith   Draw     draw = (Draw) viewer;
274cddf8d76SBarry Smith   DrawButton  button;
275cddf8d76SBarry Smith 
276416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
277416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
278416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
279416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
280cddf8d76SBarry Smith   color = DRAW_BLUE;
281416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
282cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
283416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
284cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
285cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
286cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
287cddf8d76SBarry Smith #else
288cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
289cddf8d76SBarry Smith #endif
290cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
291cddf8d76SBarry Smith     }
292cddf8d76SBarry Smith   }
293cddf8d76SBarry Smith   color = DRAW_CYAN;
294cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
295cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
296cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
297cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
298cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
299cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
300cddf8d76SBarry Smith     }
301cddf8d76SBarry Smith   }
302cddf8d76SBarry Smith   color = DRAW_RED;
303cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
304cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
305cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
306cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
307cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
308cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
309cddf8d76SBarry Smith #else
310cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
311cddf8d76SBarry Smith #endif
312cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
313416022c9SBarry Smith     }
314416022c9SBarry Smith   }
315416022c9SBarry Smith   DrawFlush(draw);
316cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
317cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
318cddf8d76SBarry Smith 
319cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
320cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
321cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
322cddf8d76SBarry Smith     DrawClear(draw);
323cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
324cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
325cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
326cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
327cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
328cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
329cddf8d76SBarry Smith     w *= scale; h *= scale;
330cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
331cddf8d76SBarry Smith     color = DRAW_BLUE;
332cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
333cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
334cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
335cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
336cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
337cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
338cddf8d76SBarry Smith #else
339cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
340cddf8d76SBarry Smith #endif
341cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
342cddf8d76SBarry Smith       }
343cddf8d76SBarry Smith     }
344cddf8d76SBarry Smith     color = DRAW_CYAN;
345cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
346cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
347cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
348cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
349cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
350cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
351cddf8d76SBarry Smith       }
352cddf8d76SBarry Smith     }
353cddf8d76SBarry Smith     color = DRAW_RED;
354cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
355cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
356cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
357cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
358cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
359cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
360cddf8d76SBarry Smith #else
361cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
362cddf8d76SBarry Smith #endif
363cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
364cddf8d76SBarry Smith       }
365cddf8d76SBarry Smith     }
366cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
367cddf8d76SBarry Smith   }
368416022c9SBarry Smith   return 0;
369416022c9SBarry Smith }
370416022c9SBarry Smith 
371416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
372416022c9SBarry Smith {
373416022c9SBarry Smith   Mat         A = (Mat) obj;
374416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
375416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
376416022c9SBarry Smith 
377416022c9SBarry Smith   if (!viewer) {
378416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
379416022c9SBarry Smith   }
380416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
381416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
382416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
383416022c9SBarry Smith     }
384416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
385416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
386416022c9SBarry Smith     }
387416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
388416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
389416022c9SBarry Smith     }
390416022c9SBarry Smith   }
391416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
392416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
393416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
39417ab2063SBarry Smith   }
39517ab2063SBarry Smith   return 0;
39617ab2063SBarry Smith }
397c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
398416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
39917ab2063SBarry Smith {
400416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
40141c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
402416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
403416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
40417ab2063SBarry Smith 
40517ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
40617ab2063SBarry Smith 
40717ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
408416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
40917ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
41017ab2063SBarry Smith     if (fshift) {
411416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
41217ab2063SBarry Smith       N = ailen[i];
41317ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
41417ab2063SBarry Smith         ip[j-fshift] = ip[j];
41517ab2063SBarry Smith         ap[j-fshift] = ap[j];
41617ab2063SBarry Smith       }
41717ab2063SBarry Smith     }
41817ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
41917ab2063SBarry Smith   }
42017ab2063SBarry Smith   if (m) {
42117ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
42217ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
42317ab2063SBarry Smith   }
42417ab2063SBarry Smith   /* reset ilen and imax for each row */
42517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
42617ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
42717ab2063SBarry Smith   }
428416022c9SBarry Smith   a->nz = ai[m] + shift;
42917ab2063SBarry Smith 
43017ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
431416022c9SBarry Smith   if (fshift && a->diag) {
4320452661fSBarry Smith     PetscFree(a->diag);
433416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
434416022c9SBarry Smith     a->diag = 0;
43517ab2063SBarry Smith   }
43676dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
43741c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
43817ab2063SBarry Smith   return 0;
43917ab2063SBarry Smith }
44017ab2063SBarry Smith 
441416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
44217ab2063SBarry Smith {
443416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
444cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
44517ab2063SBarry Smith   return 0;
44617ab2063SBarry Smith }
447416022c9SBarry Smith 
44817ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
44917ab2063SBarry Smith {
450416022c9SBarry Smith   Mat        A  = (Mat) obj;
451416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
452d5d45c9bSBarry Smith 
45317ab2063SBarry Smith #if defined(PETSC_LOG)
454416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
45517ab2063SBarry Smith #endif
4560452661fSBarry Smith   PetscFree(a->a);
4570452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4580452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
4590452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4600452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
4610452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
46276dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
4630452661fSBarry Smith   PetscFree(a);
464416022c9SBarry Smith   PLogObjectDestroy(A);
4650452661fSBarry Smith   PetscHeaderDestroy(A);
46617ab2063SBarry Smith   return 0;
46717ab2063SBarry Smith }
46817ab2063SBarry Smith 
469416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
47017ab2063SBarry Smith {
47117ab2063SBarry Smith   return 0;
47217ab2063SBarry Smith }
47317ab2063SBarry Smith 
474416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
47517ab2063SBarry Smith {
476416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
477416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
4784b0e389bSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
479416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
480416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
481416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
482e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
483e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
484e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
485e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
486df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
487df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
4884d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
4896c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
4906c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
4916c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
4926c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
4936c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
494e2f28af5SBarry Smith   else
4954d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
49617ab2063SBarry Smith   return 0;
49717ab2063SBarry Smith }
49817ab2063SBarry Smith 
499416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
50017ab2063SBarry Smith {
501416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
502416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
50317ab2063SBarry Smith   Scalar     *x, zero = 0.0;
50417ab2063SBarry Smith 
50517ab2063SBarry Smith   VecSet(&zero,v);
50617ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
507416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
508416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
509416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
510416022c9SBarry Smith       if (a->j[j]+shift == i) {
511416022c9SBarry Smith         x[i] = a->a[j];
51217ab2063SBarry Smith         break;
51317ab2063SBarry Smith       }
51417ab2063SBarry Smith     }
51517ab2063SBarry Smith   }
51617ab2063SBarry Smith   return 0;
51717ab2063SBarry Smith }
51817ab2063SBarry Smith 
51917ab2063SBarry Smith /* -------------------------------------------------------*/
52017ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
52117ab2063SBarry Smith /* -------------------------------------------------------*/
522416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
52317ab2063SBarry Smith {
524416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
52517ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
526416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
52717ab2063SBarry Smith 
52817ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
529cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
53017ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
53117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
532416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
533416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
534416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
53517ab2063SBarry Smith     alpha = x[i];
53617ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
53717ab2063SBarry Smith   }
538416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
53917ab2063SBarry Smith   return 0;
54017ab2063SBarry Smith }
541d5d45c9bSBarry Smith 
542416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
54317ab2063SBarry Smith {
544416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
54517ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
546416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
54717ab2063SBarry Smith 
54817ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
54917ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
55017ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
55117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
552416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
553416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
554416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
55517ab2063SBarry Smith     alpha = x[i];
55617ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
55717ab2063SBarry Smith   }
55817ab2063SBarry Smith   return 0;
55917ab2063SBarry Smith }
56017ab2063SBarry Smith 
561416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
56217ab2063SBarry Smith {
563416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
56417ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
565416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
56617ab2063SBarry Smith 
56717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
56817ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
569416022c9SBarry Smith   idx  = a->j;
570416022c9SBarry Smith   v    = a->a;
571416022c9SBarry Smith   ii   = a->i;
57217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
573416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
57417ab2063SBarry Smith     sum  = 0.0;
57517ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
57617ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
577416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
57817ab2063SBarry Smith     y[i] = sum;
57917ab2063SBarry Smith   }
580416022c9SBarry Smith   PLogFlops(2*a->nz - m);
58117ab2063SBarry Smith   return 0;
58217ab2063SBarry Smith }
58317ab2063SBarry Smith 
584416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
58517ab2063SBarry Smith {
586416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
58717ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
588cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
58917ab2063SBarry Smith 
59017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
59117ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
592cddf8d76SBarry Smith   idx  = a->j;
593cddf8d76SBarry Smith   v    = a->a;
594cddf8d76SBarry Smith   ii   = a->i;
59517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
596cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
59717ab2063SBarry Smith     sum  = y[i];
598cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
599cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
60017ab2063SBarry Smith     z[i] = sum;
60117ab2063SBarry Smith   }
602416022c9SBarry Smith   PLogFlops(2*a->nz);
60317ab2063SBarry Smith   return 0;
60417ab2063SBarry Smith }
60517ab2063SBarry Smith 
60617ab2063SBarry Smith /*
60717ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
60817ab2063SBarry Smith */
60917ab2063SBarry Smith 
610416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
61117ab2063SBarry Smith {
612416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
613416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
61417ab2063SBarry Smith 
6150452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
616416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
617416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
618416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
619416022c9SBarry Smith       if (a->j[j]+shift == i) {
62017ab2063SBarry Smith         diag[i] = j - shift;
62117ab2063SBarry Smith         break;
62217ab2063SBarry Smith       }
62317ab2063SBarry Smith     }
62417ab2063SBarry Smith   }
625416022c9SBarry Smith   a->diag = diag;
62617ab2063SBarry Smith   return 0;
62717ab2063SBarry Smith }
62817ab2063SBarry Smith 
629416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
63017ab2063SBarry Smith                            double fshift,int its,Vec xx)
63117ab2063SBarry Smith {
632416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
633416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
634d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
63517ab2063SBarry Smith 
63617ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
637416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
638416022c9SBarry Smith   diag = a->diag;
639416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
64017ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
64117ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
64217ab2063SBarry Smith     bs = b + shift;
64317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
644416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
645416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
646416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
647416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
64817ab2063SBarry Smith         sum  = b[i]*d/omega;
64917ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
65017ab2063SBarry Smith         x[i] = sum;
65117ab2063SBarry Smith     }
65217ab2063SBarry Smith     return 0;
65317ab2063SBarry Smith   }
65417ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
65517ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
65617ab2063SBarry Smith   }
657416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
65817ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
65917ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
66017ab2063SBarry Smith 
66117ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
66217ab2063SBarry Smith 
66317ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
66417ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
66517ab2063SBarry Smith     is the relaxation factor.
66617ab2063SBarry Smith     */
6670452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
66817ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
66917ab2063SBarry Smith 
67017ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
67117ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
672416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
673416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
674416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
675416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
67617ab2063SBarry Smith       sum  = b[i];
67717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
67817ab2063SBarry Smith       x[i] = omega*(sum/d);
67917ab2063SBarry Smith     }
68017ab2063SBarry Smith 
68117ab2063SBarry Smith     /*  t = b - (2*E - D)x */
682416022c9SBarry Smith     v = a->a;
68317ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
68417ab2063SBarry Smith 
68517ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
686416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
687416022c9SBarry Smith     diag = a->diag;
68817ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
689416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
690416022c9SBarry Smith       n    = diag[i] - a->i[i];
691416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
692416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
69317ab2063SBarry Smith       sum  = t[i];
69417ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
69517ab2063SBarry Smith       t[i] = omega*(sum/d);
69617ab2063SBarry Smith     }
69717ab2063SBarry Smith 
69817ab2063SBarry Smith     /*  x = x + t */
69917ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7000452661fSBarry Smith     PetscFree(t);
70117ab2063SBarry Smith     return 0;
70217ab2063SBarry Smith   }
70317ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
70417ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
70517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
706416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
707416022c9SBarry Smith         n    = diag[i] - a->i[i];
708416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
709416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
71017ab2063SBarry Smith         sum  = b[i];
71117ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
71217ab2063SBarry Smith         x[i] = omega*(sum/d);
71317ab2063SBarry Smith       }
71417ab2063SBarry Smith       xb = x;
71517ab2063SBarry Smith     }
71617ab2063SBarry Smith     else xb = b;
71717ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
71817ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
71917ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
720416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
72117ab2063SBarry Smith       }
72217ab2063SBarry Smith     }
72317ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
72417ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
725416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
726416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
727416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
728416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
72917ab2063SBarry Smith         sum  = xb[i];
73017ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
73117ab2063SBarry Smith         x[i] = omega*(sum/d);
73217ab2063SBarry Smith       }
73317ab2063SBarry Smith     }
73417ab2063SBarry Smith     its--;
73517ab2063SBarry Smith   }
73617ab2063SBarry Smith   while (its--) {
73717ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
73817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
739416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
740416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
741416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
742416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
74317ab2063SBarry Smith         sum  = b[i];
74417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
74517ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
74617ab2063SBarry Smith       }
74717ab2063SBarry Smith     }
74817ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
74917ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
750416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
751416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
752416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
753416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
75417ab2063SBarry Smith         sum  = b[i];
75517ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
75617ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
75717ab2063SBarry Smith       }
75817ab2063SBarry Smith     }
75917ab2063SBarry Smith   }
76017ab2063SBarry Smith   return 0;
76117ab2063SBarry Smith }
76217ab2063SBarry Smith 
763d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
76417ab2063SBarry Smith {
765416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
766416022c9SBarry Smith   *nz      = a->nz;
767416022c9SBarry Smith   *nzalloc = a->maxnz;
768416022c9SBarry Smith   *mem     = (int)A->mem;
76917ab2063SBarry Smith   return 0;
77017ab2063SBarry Smith }
77117ab2063SBarry Smith 
77217ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
77317ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
77417ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
77517ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
77617ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
77717ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
77817ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
77917ab2063SBarry Smith 
78017ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
78117ab2063SBarry Smith {
782416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
783416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
78417ab2063SBarry Smith 
78517ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
78617ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
78717ab2063SBarry Smith   if (diag) {
78817ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
789416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
790416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
791416022c9SBarry Smith         a->ilen[rows[i]] = 1;
792416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
793416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
79417ab2063SBarry Smith       }
79517ab2063SBarry Smith       else {
79617ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
79717ab2063SBarry Smith         CHKERRQ(ierr);
79817ab2063SBarry Smith       }
79917ab2063SBarry Smith     }
80017ab2063SBarry Smith   }
80117ab2063SBarry Smith   else {
80217ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
803416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
804416022c9SBarry Smith       a->ilen[rows[i]] = 0;
80517ab2063SBarry Smith     }
80617ab2063SBarry Smith   }
80717ab2063SBarry Smith   ISRestoreIndices(is,&rows);
80817ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
80917ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
81017ab2063SBarry Smith   return 0;
81117ab2063SBarry Smith }
81217ab2063SBarry Smith 
813416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
81417ab2063SBarry Smith {
815416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
816416022c9SBarry Smith   *m = a->m; *n = a->n;
81717ab2063SBarry Smith   return 0;
81817ab2063SBarry Smith }
81917ab2063SBarry Smith 
820416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
82117ab2063SBarry Smith {
822416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
823416022c9SBarry Smith   *m = 0; *n = a->m;
82417ab2063SBarry Smith   return 0;
82517ab2063SBarry Smith }
826416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
82717ab2063SBarry Smith {
828416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
829c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
83017ab2063SBarry Smith 
831416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
83217ab2063SBarry Smith 
833416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
834416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
83517ab2063SBarry Smith   if (idx) {
83617ab2063SBarry Smith     if (*nz) {
837416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
8380452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
83917ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
84017ab2063SBarry Smith     }
84117ab2063SBarry Smith     else *idx = 0;
84217ab2063SBarry Smith   }
84317ab2063SBarry Smith   return 0;
84417ab2063SBarry Smith }
84517ab2063SBarry Smith 
846416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
84717ab2063SBarry Smith {
8480452661fSBarry Smith   if (idx) {if (*idx) PetscFree(*idx);}
84917ab2063SBarry Smith   return 0;
85017ab2063SBarry Smith }
85117ab2063SBarry Smith 
852cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
85317ab2063SBarry Smith {
854416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
855416022c9SBarry Smith   Scalar     *v = a->a;
85617ab2063SBarry Smith   double     sum = 0.0;
857416022c9SBarry Smith   int        i, j,shift = a->indexshift;
85817ab2063SBarry Smith 
85917ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
860416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
86117ab2063SBarry Smith #if defined(PETSC_COMPLEX)
86217ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
86317ab2063SBarry Smith #else
86417ab2063SBarry Smith       sum += (*v)*(*v); v++;
86517ab2063SBarry Smith #endif
86617ab2063SBarry Smith     }
86717ab2063SBarry Smith     *norm = sqrt(sum);
86817ab2063SBarry Smith   }
86917ab2063SBarry Smith   else if (type == NORM_1) {
87017ab2063SBarry Smith     double *tmp;
871416022c9SBarry Smith     int    *jj = a->j;
8720452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
873cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
87417ab2063SBarry Smith     *norm = 0.0;
875416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
876a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
87717ab2063SBarry Smith     }
878416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
87917ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
88017ab2063SBarry Smith     }
8810452661fSBarry Smith     PetscFree(tmp);
88217ab2063SBarry Smith   }
88317ab2063SBarry Smith   else if (type == NORM_INFINITY) {
88417ab2063SBarry Smith     *norm = 0.0;
885416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
886416022c9SBarry Smith       v = a->a + a->i[j] + shift;
88717ab2063SBarry Smith       sum = 0.0;
888416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
889cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
89017ab2063SBarry Smith       }
89117ab2063SBarry Smith       if (sum > *norm) *norm = sum;
89217ab2063SBarry Smith     }
89317ab2063SBarry Smith   }
89417ab2063SBarry Smith   else {
89548d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
89617ab2063SBarry Smith   }
89717ab2063SBarry Smith   return 0;
89817ab2063SBarry Smith }
89917ab2063SBarry Smith 
900416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
90117ab2063SBarry Smith {
902416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
903416022c9SBarry Smith   Mat        C;
904416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
905416022c9SBarry Smith   int        shift = a->indexshift;
906d5d45c9bSBarry Smith   Scalar     *array = a->a;
90717ab2063SBarry Smith 
908416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
9090452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
910cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
91117ab2063SBarry Smith   if (shift) {
91217ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
91317ab2063SBarry Smith   }
91417ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
915416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
9160452661fSBarry Smith   PetscFree(col);
91717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
91817ab2063SBarry Smith     len = ai[i+1]-ai[i];
919416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
92017ab2063SBarry Smith     array += len; aj += len;
92117ab2063SBarry Smith   }
92217ab2063SBarry Smith   if (shift) {
92317ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
92417ab2063SBarry Smith   }
92517ab2063SBarry Smith 
926416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
927416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
92817ab2063SBarry Smith 
929416022c9SBarry Smith   if (B) {
930416022c9SBarry Smith     *B = C;
93117ab2063SBarry Smith   } else {
932416022c9SBarry Smith     /* This isn't really an in-place transpose */
9330452661fSBarry Smith     PetscFree(a->a);
9340452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
9350452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
9360452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
9370452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
9380452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
9390452661fSBarry Smith     PetscFree(a);
940416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
9410452661fSBarry Smith     PetscHeaderDestroy(C);
94217ab2063SBarry Smith   }
94317ab2063SBarry Smith   return 0;
94417ab2063SBarry Smith }
94517ab2063SBarry Smith 
946f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
94717ab2063SBarry Smith {
948416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
94917ab2063SBarry Smith   Scalar     *l,*r,x,*v;
950d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
95117ab2063SBarry Smith 
95217ab2063SBarry Smith   if (ll) {
95317ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
954f0b747eeSBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
955416022c9SBarry Smith     v = a->a;
95617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
95717ab2063SBarry Smith       x = l[i];
958416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
95917ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
96017ab2063SBarry Smith     }
96117ab2063SBarry Smith   }
96217ab2063SBarry Smith   if (rr) {
96317ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
964f0b747eeSBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
965416022c9SBarry Smith     v = a->a; jj = a->j;
96617ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
96717ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
96817ab2063SBarry Smith     }
96917ab2063SBarry Smith   }
97017ab2063SBarry Smith   return 0;
97117ab2063SBarry Smith }
97217ab2063SBarry Smith 
973cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
97417ab2063SBarry Smith {
975db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
97602834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
977a2744918SBarry Smith   register int sum,lensi;
97802834360SBarry Smith   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
979a2744918SBarry Smith   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
980db02288aSLois Curfman McInnes   Scalar       *vwork,*a_new;
981416022c9SBarry Smith   Mat          C;
98217ab2063SBarry Smith 
98317ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
98417ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
98517ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
98617ab2063SBarry Smith 
98702834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
98802834360SBarry Smith     /* special case of contiguous rows */
9890452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
99002834360SBarry Smith     starts = lens + ncols;
99102834360SBarry Smith     /* loop over new rows determining lens and starting points */
99202834360SBarry Smith     for (i=0; i<nrows; i++) {
993a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
994a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
99502834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
996d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
99702834360SBarry Smith           starts[i] = k;
99802834360SBarry Smith           break;
99902834360SBarry Smith 	}
100002834360SBarry Smith       }
1001a2744918SBarry Smith       sum = 0;
100202834360SBarry Smith       while (k < kend) {
1003d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1004a2744918SBarry Smith         sum++;
100502834360SBarry Smith       }
1006a2744918SBarry Smith       lens[i] = sum;
100702834360SBarry Smith     }
100802834360SBarry Smith     /* create submatrix */
1009cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
101008480c60SBarry Smith       int n_cols,n_rows;
101108480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
101208480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1013d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
101408480c60SBarry Smith       C = *B;
101508480c60SBarry Smith     }
101608480c60SBarry Smith     else {
101702834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
101808480c60SBarry Smith     }
1019db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1020db02288aSLois Curfman McInnes 
102102834360SBarry Smith     /* loop over rows inserting into submatrix */
1022db02288aSLois Curfman McInnes     a_new    = c->a;
1023db02288aSLois Curfman McInnes     j_new    = c->j;
1024db02288aSLois Curfman McInnes     i_new    = c->i;
1025db02288aSLois Curfman McInnes     i_new[0] = -shift;
102602834360SBarry Smith     for (i=0; i<nrows; i++) {
1027a2744918SBarry Smith       ii    = starts[i];
1028a2744918SBarry Smith       lensi = lens[i];
1029a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1030a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
103102834360SBarry Smith       }
1032a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1033a2744918SBarry Smith       a_new      += lensi;
1034a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1035a2744918SBarry Smith       c->ilen[i]  = lensi;
103602834360SBarry Smith     }
10370452661fSBarry Smith     PetscFree(lens);
103802834360SBarry Smith   }
103902834360SBarry Smith   else {
104002834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10410452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
104202834360SBarry Smith     ssmap = smap + shift;
10437eb43aa7SLois Curfman McInnes     cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork);
104402834360SBarry Smith     lens  = cwork + ncols;
10450452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1046cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
104717ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
104802834360SBarry Smith     /* determine lens of each row */
104902834360SBarry Smith     for (i=0; i<nrows; i++) {
1050d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
105102834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
105202834360SBarry Smith       lens[i] = 0;
105302834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1054d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
105502834360SBarry Smith           lens[i]++;
105602834360SBarry Smith         }
105702834360SBarry Smith       }
105802834360SBarry Smith     }
105917ab2063SBarry Smith     /* Create and fill new matrix */
1060a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
106108480c60SBarry Smith       int n_cols,n_rows;
106208480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
106308480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1064d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
106508480c60SBarry Smith       C = *B;
106608480c60SBarry Smith     }
106708480c60SBarry Smith     else {
106802834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
106908480c60SBarry Smith     }
107017ab2063SBarry Smith     for (i=0; i<nrows; i++) {
107117ab2063SBarry Smith       nznew  = 0;
1072d8ced48eSBarry Smith       kstart = ai[irow[i]]+shift;
1073416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
107417ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
107502834360SBarry Smith         if (ssmap[a->j[k]]) {
107602834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1077416022c9SBarry Smith           vwork[nznew++] = a->a[k];
107817ab2063SBarry Smith         }
107917ab2063SBarry Smith       }
1080416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
108117ab2063SBarry Smith     }
108202834360SBarry Smith     /* Free work space */
108302834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10840452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
108502834360SBarry Smith   }
1086416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1087416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
108817ab2063SBarry Smith 
108917ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1090416022c9SBarry Smith   *B = C;
109117ab2063SBarry Smith   return 0;
109217ab2063SBarry Smith }
109317ab2063SBarry Smith 
1094a871dcd8SBarry Smith /*
109563b91edcSBarry Smith      note: This can only work for identity for row and col. It would
109663b91edcSBarry Smith    be good to check this and otherwise generate an error.
1097a871dcd8SBarry Smith */
109863b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1099a871dcd8SBarry Smith {
110063b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
110108480c60SBarry Smith   int        ierr;
110263b91edcSBarry Smith   Mat        outA;
110363b91edcSBarry Smith 
1104a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1105a871dcd8SBarry Smith 
110663b91edcSBarry Smith   outA          = inA;
110763b91edcSBarry Smith   inA->factor   = FACTOR_LU;
110863b91edcSBarry Smith   a->row        = row;
110963b91edcSBarry Smith   a->col        = col;
111063b91edcSBarry Smith 
11110452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
111263b91edcSBarry Smith 
111308480c60SBarry Smith   if (!a->diag) {
111408480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
111563b91edcSBarry Smith   }
111663b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1117a871dcd8SBarry Smith   return 0;
1118a871dcd8SBarry Smith }
1119a871dcd8SBarry Smith 
1120f0b747eeSBarry Smith #include "pinclude/plapack.h"
1121f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1122f0b747eeSBarry Smith {
1123f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1124f0b747eeSBarry Smith   int        one = 1;
1125f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1126f0b747eeSBarry Smith   PLogFlops(a->nz);
1127f0b747eeSBarry Smith   return 0;
1128f0b747eeSBarry Smith }
1129f0b747eeSBarry Smith 
1130cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1131cddf8d76SBarry Smith                                     Mat **B)
1132cddf8d76SBarry Smith {
1133cddf8d76SBarry Smith   int ierr,i;
1134cddf8d76SBarry Smith 
1135cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
11360452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1137cddf8d76SBarry Smith   }
1138cddf8d76SBarry Smith 
1139cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1140cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1141cddf8d76SBarry Smith   }
1142cddf8d76SBarry Smith   return 0;
1143cddf8d76SBarry Smith }
1144cddf8d76SBarry Smith 
1145e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
11464dcbc457SBarry Smith {
1147e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11488a047759SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *table, *nidx, isz, val;
11498a047759SSatish Balay   int        start, end, *ai, *aj;
11504dcbc457SBarry Smith 
11518a047759SSatish Balay   shift = a->indexshift;
1152e4d965acSSatish Balay   m     = a->m;
1153e4d965acSSatish Balay   ai    = a->i;
11548a047759SSatish Balay   aj    = a->j+shift;
11558a047759SSatish Balay 
115604a348a9SBarry Smith   table = (int *) PetscMalloc(2*m*sizeof(int)); CHKPTRQ(table); nidx = table + m;
11578a047759SSatish Balay 
1158e4d965acSSatish Balay   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
1159e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1160e4d965acSSatish Balay     /* Initialise the two local arrays */
1161e4d965acSSatish Balay     isz  = 0;
1162e4d965acSSatish Balay     PetscMemzero(table,m*sizeof(int));
1163e4d965acSSatish Balay 
1164e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
11654dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
1166e4d965acSSatish Balay     ierr = ISGetLocalSize(is[i],&n);  CHKERRQ(ierr);
1167e4d965acSSatish Balay 
1168e4d965acSSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1169e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
1170e4d965acSSatish Balay       if(!table[idx[j]]++) { nidx[isz++] = idx[j];}
11714dcbc457SBarry Smith     }
1172e4d965acSSatish Balay 
117304a348a9SBarry Smith     k = 0;
117404a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
117504a348a9SBarry Smith       n = isz;
117604a348a9SBarry Smith       for ( ; k<n ; k++){ /* do only those rows in nidx[k] not yet done */
1177e4d965acSSatish Balay         row   = nidx[k];
1178e4d965acSSatish Balay         start = ai[row];
1179e4d965acSSatish Balay         end   = ai[row+1];
118004a348a9SBarry Smith         for ( l = start; l<end ; l++){
11818a047759SSatish Balay           val = aj[l] + shift;
1182e4d965acSSatish Balay           if (!table[val]++) {nidx[isz++] = val;}
1183e4d965acSSatish Balay         }
1184e4d965acSSatish Balay       }
1185e4d965acSSatish Balay     }
1186e4d965acSSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
1187e4d965acSSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1188e4d965acSSatish Balay     ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1189e4d965acSSatish Balay   }
119004a348a9SBarry Smith   PetscFree(table);
1191e4d965acSSatish Balay   return 0;
11924dcbc457SBarry Smith }
119317ab2063SBarry Smith 
1194682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1195682d7d0cSBarry Smith {
1196682d7d0cSBarry Smith   static int called = 0;
1197682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1198682d7d0cSBarry Smith 
1199682d7d0cSBarry Smith   if (called) return 0; else called = 1;
12002a7368beSLois Curfman McInnes   MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1201682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
12022a7368beSLois Curfman McInnes   MPIU_printf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
1203682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
1204682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1205682d7d0cSBarry Smith #if defined(HAVE_ESSL)
1206682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1207682d7d0cSBarry Smith #endif
1208682d7d0cSBarry Smith   return 0;
1209682d7d0cSBarry Smith }
1210682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
121117ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
121217ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1213416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1214416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
121517ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
121617ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
121717ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
121817ab2063SBarry Smith        MatRelax_SeqAIJ,
121917ab2063SBarry Smith        MatTranspose_SeqAIJ,
122017ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
1221f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
122217ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
122317ab2063SBarry Smith        MatCompress_SeqAIJ,
122417ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
122517ab2063SBarry Smith        MatGetReordering_SeqAIJ,
122617ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
122717ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
122817ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
122917ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1230416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
12313d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1232cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
12337eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1234682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1235f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
1236f0b747eeSBarry Smith        MatScale_SeqAIJ};
123717ab2063SBarry Smith 
123817ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
123917ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
124017ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
124117ab2063SBarry Smith 
124217ab2063SBarry Smith /*@C
1243682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
124417ab2063SBarry Smith                      (the default uniprocessor PETSc format).
124517ab2063SBarry Smith 
124617ab2063SBarry Smith    Input Parameters:
124717ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
124817ab2063SBarry Smith .  m - number of rows
124917ab2063SBarry Smith .  n - number of columns
125017ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
125117ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
125217ab2063SBarry Smith 
125317ab2063SBarry Smith    Output Parameter:
1254416022c9SBarry Smith .  A - the matrix
125517ab2063SBarry Smith 
125617ab2063SBarry Smith    Notes:
125717ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
125817ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
12590002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
12600002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
126117ab2063SBarry Smith 
126217ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1263a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
126417ab2063SBarry Smith    allocation.
126517ab2063SBarry Smith 
1266682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1267682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1268682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
12696c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
12706c7ebb05SLois Curfman McInnes 
12716c7ebb05SLois Curfman McInnes    Options Database Keys:
12726c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
12736c7ebb05SLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)
12744b0e389bSBarry Smith $    -mat_aij_oneindex - Internally using indexing starting at 1 rather then zero.
12754b0e389bSBarry Smith .                        Note: You still index entries starting at 0!
127617ab2063SBarry Smith 
127717ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
127817ab2063SBarry Smith @*/
1279416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
128017ab2063SBarry Smith {
1281416022c9SBarry Smith   Mat        B;
1282416022c9SBarry Smith   Mat_SeqAIJ *b;
128369957df2SSatish Balay   int        i,len,ierr, flg;
1284d5d45c9bSBarry Smith 
1285416022c9SBarry Smith   *A      = 0;
12860452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1287416022c9SBarry Smith   PLogObjectCreate(B);
12880452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1289416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1290416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1291416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1292416022c9SBarry Smith   B->factor           = 0;
1293416022c9SBarry Smith   B->lupivotthreshold = 1.0;
129469957df2SSatish Balay   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, \
129569957df2SSatish Balay                           &flg); CHKERRQ(ierr);
1296416022c9SBarry Smith   b->row              = 0;
1297416022c9SBarry Smith   b->col              = 0;
1298416022c9SBarry Smith   b->indexshift       = 0;
129969957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
130069957df2SSatish Balay   if (flg) b->indexshift = -1;
130117ab2063SBarry Smith 
1302416022c9SBarry Smith   b->m       = m;
1303416022c9SBarry Smith   b->n       = n;
13040452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1305b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
13067b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
13077b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1308416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
130917ab2063SBarry Smith     nz = nz*m;
131017ab2063SBarry Smith   }
131117ab2063SBarry Smith   else {
131217ab2063SBarry Smith     nz = 0;
1313416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
131417ab2063SBarry Smith   }
131517ab2063SBarry Smith 
131617ab2063SBarry Smith   /* allocate the matrix space */
1317416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
13180452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1319416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1320cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1321416022c9SBarry Smith   b->i  = b->j + nz;
1322416022c9SBarry Smith   b->singlemalloc = 1;
132317ab2063SBarry Smith 
1324416022c9SBarry Smith   b->i[0] = -b->indexshift;
132517ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1326416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
132717ab2063SBarry Smith   }
132817ab2063SBarry Smith 
1329416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
13300452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1331416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1332416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
133317ab2063SBarry Smith 
1334416022c9SBarry Smith   b->nz               = 0;
1335416022c9SBarry Smith   b->maxnz            = nz;
1336416022c9SBarry Smith   b->sorted           = 0;
1337416022c9SBarry Smith   b->roworiented      = 1;
1338416022c9SBarry Smith   b->nonew            = 0;
1339416022c9SBarry Smith   b->diag             = 0;
1340416022c9SBarry Smith   b->solve_work       = 0;
134171bd300dSLois Curfman McInnes   b->spptr            = 0;
1342754ec7b1SSatish Balay   b->inode.node_count = 0;
1343754ec7b1SSatish Balay   b->inode.size       = 0;
13446c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
13456c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
134617ab2063SBarry Smith 
1347416022c9SBarry Smith   *A = B;
1348*4b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
1349*4b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
135069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
135169957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
1352*4b14c69eSBarry Smith #endif
135369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
135469957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
135569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
135669957df2SSatish Balay   if (flg) {
1357416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1358416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
135917ab2063SBarry Smith   }
136069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
136169957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
136217ab2063SBarry Smith   return 0;
136317ab2063SBarry Smith }
136417ab2063SBarry Smith 
13653d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
136617ab2063SBarry Smith {
1367416022c9SBarry Smith   Mat        C;
1368416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
136908480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
137017ab2063SBarry Smith 
13714043dd9cSLois Curfman McInnes   *B = 0;
13720452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1373416022c9SBarry Smith   PLogObjectCreate(C);
13740452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
137541c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1376416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1377416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1378416022c9SBarry Smith   C->factor     = A->factor;
1379416022c9SBarry Smith   c->row        = 0;
1380416022c9SBarry Smith   c->col        = 0;
1381416022c9SBarry Smith   c->indexshift = shift;
1382c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
138317ab2063SBarry Smith 
1384416022c9SBarry Smith   c->m          = a->m;
1385416022c9SBarry Smith   c->n          = a->n;
138617ab2063SBarry Smith 
13870452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
13880452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
138917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1390416022c9SBarry Smith     c->imax[i] = a->imax[i];
1391416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
139217ab2063SBarry Smith   }
139317ab2063SBarry Smith 
139417ab2063SBarry Smith   /* allocate the matrix space */
1395416022c9SBarry Smith   c->singlemalloc = 1;
1396416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
13970452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1398416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1399416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1400416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
140117ab2063SBarry Smith   if (m > 0) {
1402416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
140308480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1404416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
140517ab2063SBarry Smith     }
140608480c60SBarry Smith   }
140717ab2063SBarry Smith 
1408416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1409416022c9SBarry Smith   c->sorted      = a->sorted;
1410416022c9SBarry Smith   c->roworiented = a->roworiented;
1411416022c9SBarry Smith   c->nonew       = a->nonew;
141217ab2063SBarry Smith 
1413416022c9SBarry Smith   if (a->diag) {
14140452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1415416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
141617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1417416022c9SBarry Smith       c->diag[i] = a->diag[i];
141817ab2063SBarry Smith     }
141917ab2063SBarry Smith   }
1420416022c9SBarry Smith   else c->diag          = 0;
14216c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
14226c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1423754ec7b1SSatish Balay   if (a->inode.size){
1424754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1425754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1426754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1427754ec7b1SSatish Balay   } else {
1428754ec7b1SSatish Balay     c->inode.size       = 0;
1429754ec7b1SSatish Balay     c->inode.node_count = 0;
1430754ec7b1SSatish Balay   }
1431416022c9SBarry Smith   c->nz                 = a->nz;
1432416022c9SBarry Smith   c->maxnz              = a->maxnz;
1433416022c9SBarry Smith   c->solve_work         = 0;
143476dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1435754ec7b1SSatish Balay 
1436416022c9SBarry Smith   *B = C;
143717ab2063SBarry Smith   return 0;
143817ab2063SBarry Smith }
143917ab2063SBarry Smith 
1440416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
144117ab2063SBarry Smith {
1442416022c9SBarry Smith   Mat_SeqAIJ   *a;
1443416022c9SBarry Smith   Mat          B;
144417699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
144517ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
144617ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
144717ab2063SBarry Smith 
144817699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
144917699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
145017ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1451416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
145248d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
145317ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
145417ab2063SBarry Smith 
145517ab2063SBarry Smith   /* read in row lengths */
14560452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1457416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
145817ab2063SBarry Smith 
145917ab2063SBarry Smith   /* create our matrix */
1460416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1461416022c9SBarry Smith   B = *A;
1462416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1463416022c9SBarry Smith   shift = a->indexshift;
146417ab2063SBarry Smith 
146517ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1466416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
146717ab2063SBarry Smith   if (shift) {
146817ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1469416022c9SBarry Smith       a->j[i] += 1;
147017ab2063SBarry Smith     }
147117ab2063SBarry Smith   }
147217ab2063SBarry Smith 
147317ab2063SBarry Smith   /* read in nonzero values */
1474416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
147517ab2063SBarry Smith 
147617ab2063SBarry Smith   /* set matrix "i" values */
1477416022c9SBarry Smith   a->i[0] = -shift;
147817ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1479416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1480416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
148117ab2063SBarry Smith   }
14820452661fSBarry Smith   PetscFree(rowlengths);
148317ab2063SBarry Smith 
1484416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1485416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
148617ab2063SBarry Smith   return 0;
148717ab2063SBarry Smith }
148817ab2063SBarry Smith 
148917ab2063SBarry Smith 
149017ab2063SBarry Smith 
1491