xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 3d1612f7f6a78b84146aff229cce396c8cf5cec5)
117ab2063SBarry Smith #ifndef lint
2*3d1612f7SBarry Smith static char vcid[] = "$Id: aij.c,v 1.126 1995/12/21 18:31:37 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"
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;
64416022c9SBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen;
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");
78416022c9SBarry Smith       col = in[l] - shift; value = *v++;
79416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
80416022c9SBarry Smith       while (high-low > 5) {
81416022c9SBarry Smith         t = (low+high)/2;
82416022c9SBarry Smith         if (rp[t] > col) high = t;
83416022c9SBarry Smith         else             low  = t;
8417ab2063SBarry Smith       }
85416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
8617ab2063SBarry Smith         if (rp[i] > col) break;
8717ab2063SBarry Smith         if (rp[i] == col) {
88416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
8917ab2063SBarry Smith           else                  ap[i] = value;
9017ab2063SBarry Smith           goto noinsert;
9117ab2063SBarry Smith         }
9217ab2063SBarry Smith       }
9317ab2063SBarry Smith       if (nonew) goto noinsert;
9417ab2063SBarry Smith       if (nrow >= rmax) {
9517ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
96416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
9717ab2063SBarry Smith         Scalar *new_a;
9817ab2063SBarry Smith 
9917ab2063SBarry Smith         /* malloc new storage space */
100416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1010452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
10217ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
10317ab2063SBarry Smith         new_i   = new_j + new_nz;
10417ab2063SBarry Smith 
10517ab2063SBarry Smith         /* copy over old data into new slots */
10617ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
107416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
108416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
109416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
110416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
11117ab2063SBarry Smith                                                            len*sizeof(int));
112416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
113416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
11417ab2063SBarry Smith                                                            len*sizeof(Scalar));
11517ab2063SBarry Smith         /* free up old matrix storage */
1160452661fSBarry Smith         PetscFree(a->a);
1170452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
118416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
119416022c9SBarry Smith         a->singlemalloc = 1;
12017ab2063SBarry Smith 
12117ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
122416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
123416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
124416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
12517ab2063SBarry Smith       }
126416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
127416022c9SBarry Smith       /* shift up all the later entries in this row */
128416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
12917ab2063SBarry Smith         rp[ii+1] = rp[ii];
13017ab2063SBarry Smith         ap[ii+1] = ap[ii];
13117ab2063SBarry Smith       }
13217ab2063SBarry Smith       rp[i] = col;
13317ab2063SBarry Smith       ap[i] = value;
13417ab2063SBarry Smith       noinsert:;
135416022c9SBarry Smith       low = i + 1;
13617ab2063SBarry Smith     }
13717ab2063SBarry Smith     ailen[row] = nrow;
13817ab2063SBarry Smith   }
13917ab2063SBarry Smith   return 0;
14017ab2063SBarry Smith }
14117ab2063SBarry Smith 
1427eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1437eb43aa7SLois Curfman McInnes {
1447eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
145b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1467eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1477eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
1487eb43aa7SLois Curfman McInnes 
1497eb43aa7SLois Curfman McInnes   if (!a->assembled) SETERRQ(1,"MatGetValues_SeqAIJ:Not for unassembled matrix");
1507eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
1517eb43aa7SLois Curfman McInnes     row  = im[k];
1527eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
1537eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
1547eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
1557eb43aa7SLois Curfman McInnes     nrow = ailen[row];
1567eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
1577eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
1587eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
1597eb43aa7SLois Curfman McInnes       col = in[l] - shift;
1607eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
1617eb43aa7SLois Curfman McInnes       while (high-low > 5) {
1627eb43aa7SLois Curfman McInnes         t = (low+high)/2;
1637eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
1647eb43aa7SLois Curfman McInnes         else             low  = t;
1657eb43aa7SLois Curfman McInnes       }
1667eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
1677eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
1687eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
169b49de8d1SLois Curfman McInnes           *v++ = ap[i];
1707eb43aa7SLois Curfman McInnes           goto finished;
1717eb43aa7SLois Curfman McInnes         }
1727eb43aa7SLois Curfman McInnes       }
173b49de8d1SLois Curfman McInnes       *v++ = zero;
1747eb43aa7SLois Curfman McInnes       finished:;
1757eb43aa7SLois Curfman McInnes     }
1767eb43aa7SLois Curfman McInnes   }
1777eb43aa7SLois Curfman McInnes   return 0;
1787eb43aa7SLois Curfman McInnes }
1797eb43aa7SLois Curfman McInnes 
18017ab2063SBarry Smith #include "draw.h"
18117ab2063SBarry Smith #include "pinclude/pviewer.h"
182416022c9SBarry Smith #include "sysio.h"
18317ab2063SBarry Smith 
184416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
18517ab2063SBarry Smith {
186416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
187416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
18817ab2063SBarry Smith 
189416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1900452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
191416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
192416022c9SBarry Smith   col_lens[1] = a->m;
193416022c9SBarry Smith   col_lens[2] = a->n;
194416022c9SBarry Smith   col_lens[3] = a->nz;
195416022c9SBarry Smith 
196416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
197416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
198416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
19917ab2063SBarry Smith   }
200416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
2010452661fSBarry Smith   PetscFree(col_lens);
202416022c9SBarry Smith 
203416022c9SBarry Smith   /* store column indices (zero start index) */
204416022c9SBarry Smith   if (a->indexshift) {
205416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
20617ab2063SBarry Smith   }
207416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
208416022c9SBarry Smith   if (a->indexshift) {
209416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
21017ab2063SBarry Smith   }
211416022c9SBarry Smith 
212416022c9SBarry Smith   /* store nonzero values */
213416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
21417ab2063SBarry Smith   return 0;
21517ab2063SBarry Smith }
216416022c9SBarry Smith 
217416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
218416022c9SBarry Smith {
219416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
220416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
22117ab2063SBarry Smith   FILE        *fd;
22217ab2063SBarry Smith   char        *outputname;
22317ab2063SBarry Smith 
224416022c9SBarry Smith   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
225416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
226416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
22717ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
22808480c60SBarry Smith     /* no need to print additional information */ ;
22917ab2063SBarry Smith   }
23017ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
23117ab2063SBarry Smith     int nz, nzalloc, mem;
232416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
233416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
23417ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
23517ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
23617ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
23717ab2063SBarry Smith 
23817ab2063SBarry Smith     for (i=0; i<m; i++) {
239416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
24017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
241416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
242416022c9SBarry Smith                    imag(a->a[j]));
24317ab2063SBarry Smith #else
244416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
24517ab2063SBarry Smith #endif
24617ab2063SBarry Smith       }
24717ab2063SBarry Smith     }
24817ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
24917ab2063SBarry Smith   }
25017ab2063SBarry Smith   else {
25117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
25217ab2063SBarry Smith       fprintf(fd,"row %d:",i);
253416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
25417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
255416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
256416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
25717ab2063SBarry Smith         }
25817ab2063SBarry Smith         else {
259416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
26017ab2063SBarry Smith         }
26117ab2063SBarry Smith #else
262416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
26317ab2063SBarry Smith #endif
26417ab2063SBarry Smith       }
26517ab2063SBarry Smith       fprintf(fd,"\n");
26617ab2063SBarry Smith     }
26717ab2063SBarry Smith   }
26817ab2063SBarry Smith   fflush(fd);
269416022c9SBarry Smith   return 0;
270416022c9SBarry Smith }
271416022c9SBarry Smith 
272416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
273416022c9SBarry Smith {
274416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
275cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
276cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
277d7e8b826SBarry Smith   Draw     draw = (Draw) viewer;
278cddf8d76SBarry Smith   DrawButton  button;
279cddf8d76SBarry Smith 
280416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
281416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
282416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
283416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
284cddf8d76SBarry Smith   color = DRAW_BLUE;
285416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
286cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
287416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
288cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
289cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
290cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
291cddf8d76SBarry Smith #else
292cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
293cddf8d76SBarry Smith #endif
294cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
295cddf8d76SBarry Smith     }
296cddf8d76SBarry Smith   }
297cddf8d76SBarry Smith   color = DRAW_CYAN;
298cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
299cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
300cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
301cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
302cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
303cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
304cddf8d76SBarry Smith     }
305cddf8d76SBarry Smith   }
306cddf8d76SBarry Smith   color = DRAW_RED;
307cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
308cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
309cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
310cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
311cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
312cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
313cddf8d76SBarry Smith #else
314cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
315cddf8d76SBarry Smith #endif
316cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
317416022c9SBarry Smith     }
318416022c9SBarry Smith   }
319416022c9SBarry Smith   DrawFlush(draw);
320cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
321cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
322cddf8d76SBarry Smith 
323cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
324cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
325cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
326cddf8d76SBarry Smith     DrawClear(draw);
327cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
328cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
329cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
330cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
331cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
332cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
333cddf8d76SBarry Smith     w *= scale; h *= scale;
334cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
335cddf8d76SBarry Smith     color = DRAW_BLUE;
336cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
337cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
338cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
339cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
340cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
341cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
342cddf8d76SBarry Smith #else
343cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
344cddf8d76SBarry Smith #endif
345cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
346cddf8d76SBarry Smith       }
347cddf8d76SBarry Smith     }
348cddf8d76SBarry Smith     color = DRAW_CYAN;
349cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
350cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
351cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
352cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
353cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
354cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
355cddf8d76SBarry Smith       }
356cddf8d76SBarry Smith     }
357cddf8d76SBarry Smith     color = DRAW_RED;
358cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
359cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
360cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
361cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
362cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
363cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
364cddf8d76SBarry Smith #else
365cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
366cddf8d76SBarry Smith #endif
367cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
368cddf8d76SBarry Smith       }
369cddf8d76SBarry Smith     }
370cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
371cddf8d76SBarry Smith   }
372416022c9SBarry Smith   return 0;
373416022c9SBarry Smith }
374416022c9SBarry Smith 
375416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
376416022c9SBarry Smith {
377416022c9SBarry Smith   Mat         A = (Mat) obj;
378416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
379416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
380416022c9SBarry Smith 
381416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
382416022c9SBarry Smith   if (!viewer) {
383416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
384416022c9SBarry Smith   }
385416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
386416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
387416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
388416022c9SBarry Smith     }
389416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
390416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
391416022c9SBarry Smith     }
392416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
393416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
394416022c9SBarry Smith     }
395416022c9SBarry Smith   }
396416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
397416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
398416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
39917ab2063SBarry Smith   }
40017ab2063SBarry Smith   return 0;
40117ab2063SBarry Smith }
40241c01911SSatish Balay int Mat_AIJ_CheckInode(Mat);
403416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
40417ab2063SBarry Smith {
405416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
40641c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
407416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
408416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
40917ab2063SBarry Smith 
41017ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
41117ab2063SBarry Smith 
41217ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
413416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
41417ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
41517ab2063SBarry Smith     if (fshift) {
416416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
41717ab2063SBarry Smith       N = ailen[i];
41817ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
41917ab2063SBarry Smith         ip[j-fshift] = ip[j];
42017ab2063SBarry Smith         ap[j-fshift] = ap[j];
42117ab2063SBarry Smith       }
42217ab2063SBarry Smith     }
42317ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
42417ab2063SBarry Smith   }
42517ab2063SBarry Smith   if (m) {
42617ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
42717ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
42817ab2063SBarry Smith   }
42917ab2063SBarry Smith   /* reset ilen and imax for each row */
43017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
43117ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
43217ab2063SBarry Smith   }
433416022c9SBarry Smith   a->nz = ai[m] + shift;
43417ab2063SBarry Smith 
43517ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
436416022c9SBarry Smith   if (fshift && a->diag) {
4370452661fSBarry Smith     PetscFree(a->diag);
438416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
439416022c9SBarry Smith     a->diag = 0;
44017ab2063SBarry Smith   }
44176dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
44241c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
443416022c9SBarry Smith   a->assembled = 1;
44417ab2063SBarry Smith   return 0;
44517ab2063SBarry Smith }
44617ab2063SBarry Smith 
447416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
44817ab2063SBarry Smith {
449416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
450cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
45117ab2063SBarry Smith   return 0;
45217ab2063SBarry Smith }
453416022c9SBarry Smith 
45417ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
45517ab2063SBarry Smith {
456416022c9SBarry Smith   Mat        A  = (Mat) obj;
457416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
458d5d45c9bSBarry Smith 
45917ab2063SBarry Smith #if defined(PETSC_LOG)
460416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
46117ab2063SBarry Smith #endif
4620452661fSBarry Smith   PetscFree(a->a);
4630452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4640452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
4650452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4660452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
4670452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
46876dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
4690452661fSBarry Smith   PetscFree(a);
470416022c9SBarry Smith   PLogObjectDestroy(A);
4710452661fSBarry Smith   PetscHeaderDestroy(A);
47217ab2063SBarry Smith   return 0;
47317ab2063SBarry Smith }
47417ab2063SBarry Smith 
475416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
47617ab2063SBarry Smith {
47717ab2063SBarry Smith   return 0;
47817ab2063SBarry Smith }
47917ab2063SBarry Smith 
480416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
48117ab2063SBarry Smith {
482416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
483416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
484416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
485416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
486416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
487e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
488e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
489e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
490e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
491df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
492df719601SLois Curfman McInnes   else if (op == COLUMN_ORIENTED)
4934d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");}
494df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
4954d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
4966c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
4976c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
4986c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
4996c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
5006c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
501e2f28af5SBarry Smith   else
5024d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
50317ab2063SBarry Smith   return 0;
50417ab2063SBarry Smith }
50517ab2063SBarry Smith 
506416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
50717ab2063SBarry Smith {
508416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
509416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
51017ab2063SBarry Smith   Scalar     *x, zero = 0.0;
51117ab2063SBarry Smith 
512416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
51317ab2063SBarry Smith   VecSet(&zero,v);
51417ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
515416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
516416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
517416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
518416022c9SBarry Smith       if (a->j[j]+shift == i) {
519416022c9SBarry Smith         x[i] = a->a[j];
52017ab2063SBarry Smith         break;
52117ab2063SBarry Smith       }
52217ab2063SBarry Smith     }
52317ab2063SBarry Smith   }
52417ab2063SBarry Smith   return 0;
52517ab2063SBarry Smith }
52617ab2063SBarry Smith 
52717ab2063SBarry Smith /* -------------------------------------------------------*/
52817ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
52917ab2063SBarry Smith /* -------------------------------------------------------*/
530416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
53117ab2063SBarry Smith {
532416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
53317ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
534416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
53517ab2063SBarry Smith 
536416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
53717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
538cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
53917ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
54017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
541416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
542416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
543416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
54417ab2063SBarry Smith     alpha = x[i];
54517ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
54617ab2063SBarry Smith   }
547416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
54817ab2063SBarry Smith   return 0;
54917ab2063SBarry Smith }
550d5d45c9bSBarry Smith 
551416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
55217ab2063SBarry Smith {
553416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
55417ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
555416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
55617ab2063SBarry Smith 
557416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
55817ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
55917ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
56017ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
56117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
562416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
563416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
564416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
56517ab2063SBarry Smith     alpha = x[i];
56617ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
56717ab2063SBarry Smith   }
56817ab2063SBarry Smith   return 0;
56917ab2063SBarry Smith }
57017ab2063SBarry Smith 
571416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
57217ab2063SBarry Smith {
573416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
57417ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
575416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
57617ab2063SBarry Smith 
577416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
57817ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
57917ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
580416022c9SBarry Smith   idx  = a->j;
581416022c9SBarry Smith   v    = a->a;
582416022c9SBarry Smith   ii   = a->i;
58317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
584416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
58517ab2063SBarry Smith     sum  = 0.0;
58617ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
58717ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
588416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
58917ab2063SBarry Smith     y[i] = sum;
59017ab2063SBarry Smith   }
591416022c9SBarry Smith   PLogFlops(2*a->nz - m);
59217ab2063SBarry Smith   return 0;
59317ab2063SBarry Smith }
59417ab2063SBarry Smith 
595416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
59617ab2063SBarry Smith {
597416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
59817ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
599cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
60017ab2063SBarry Smith 
60148d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
60217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
60317ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
604cddf8d76SBarry Smith   idx  = a->j;
605cddf8d76SBarry Smith   v    = a->a;
606cddf8d76SBarry Smith   ii   = a->i;
60717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
608cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
60917ab2063SBarry Smith     sum  = y[i];
610cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
611cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
61217ab2063SBarry Smith     z[i] = sum;
61317ab2063SBarry Smith   }
614416022c9SBarry Smith   PLogFlops(2*a->nz);
61517ab2063SBarry Smith   return 0;
61617ab2063SBarry Smith }
61717ab2063SBarry Smith 
61817ab2063SBarry Smith /*
61917ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
62017ab2063SBarry Smith */
62117ab2063SBarry Smith 
622416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
62317ab2063SBarry Smith {
624416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
625416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
62617ab2063SBarry Smith 
627416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
6280452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
629416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
630416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
631416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
632416022c9SBarry Smith       if (a->j[j]+shift == i) {
63317ab2063SBarry Smith         diag[i] = j - shift;
63417ab2063SBarry Smith         break;
63517ab2063SBarry Smith       }
63617ab2063SBarry Smith     }
63717ab2063SBarry Smith   }
638416022c9SBarry Smith   a->diag = diag;
63917ab2063SBarry Smith   return 0;
64017ab2063SBarry Smith }
64117ab2063SBarry Smith 
642416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
64317ab2063SBarry Smith                            double fshift,int its,Vec xx)
64417ab2063SBarry Smith {
645416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
646416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
647d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
64817ab2063SBarry Smith 
64917ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
650416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
651416022c9SBarry Smith   diag = a->diag;
652416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
65317ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
65417ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
65517ab2063SBarry Smith     bs = b + shift;
65617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
657416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
658416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
659416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
660416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
66117ab2063SBarry Smith         sum  = b[i]*d/omega;
66217ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
66317ab2063SBarry Smith         x[i] = sum;
66417ab2063SBarry Smith     }
66517ab2063SBarry Smith     return 0;
66617ab2063SBarry Smith   }
66717ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
66817ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
66917ab2063SBarry Smith   }
670416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
67117ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
67217ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
67317ab2063SBarry Smith 
67417ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
67517ab2063SBarry Smith 
67617ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
67717ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
67817ab2063SBarry Smith     is the relaxation factor.
67917ab2063SBarry Smith     */
6800452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
68117ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
68217ab2063SBarry Smith 
68317ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
68417ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
685416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
686416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
687416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
688416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
68917ab2063SBarry Smith       sum  = b[i];
69017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
69117ab2063SBarry Smith       x[i] = omega*(sum/d);
69217ab2063SBarry Smith     }
69317ab2063SBarry Smith 
69417ab2063SBarry Smith     /*  t = b - (2*E - D)x */
695416022c9SBarry Smith     v = a->a;
69617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
69717ab2063SBarry Smith 
69817ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
699416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
700416022c9SBarry Smith     diag = a->diag;
70117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
702416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
703416022c9SBarry Smith       n    = diag[i] - a->i[i];
704416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
705416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
70617ab2063SBarry Smith       sum  = t[i];
70717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
70817ab2063SBarry Smith       t[i] = omega*(sum/d);
70917ab2063SBarry Smith     }
71017ab2063SBarry Smith 
71117ab2063SBarry Smith     /*  x = x + t */
71217ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7130452661fSBarry Smith     PetscFree(t);
71417ab2063SBarry Smith     return 0;
71517ab2063SBarry Smith   }
71617ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
71717ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
71817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
719416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
720416022c9SBarry Smith         n    = diag[i] - a->i[i];
721416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
722416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
72317ab2063SBarry Smith         sum  = b[i];
72417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
72517ab2063SBarry Smith         x[i] = omega*(sum/d);
72617ab2063SBarry Smith       }
72717ab2063SBarry Smith       xb = x;
72817ab2063SBarry Smith     }
72917ab2063SBarry Smith     else xb = b;
73017ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
73117ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
73217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
733416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
73417ab2063SBarry Smith       }
73517ab2063SBarry Smith     }
73617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
73717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
738416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
739416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
740416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
741416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
74217ab2063SBarry Smith         sum  = xb[i];
74317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
74417ab2063SBarry Smith         x[i] = omega*(sum/d);
74517ab2063SBarry Smith       }
74617ab2063SBarry Smith     }
74717ab2063SBarry Smith     its--;
74817ab2063SBarry Smith   }
74917ab2063SBarry Smith   while (its--) {
75017ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
75117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
752416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
753416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
754416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
755416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
75617ab2063SBarry Smith         sum  = b[i];
75717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
75817ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
75917ab2063SBarry Smith       }
76017ab2063SBarry Smith     }
76117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
76217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
763416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
764416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
765416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
766416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
76717ab2063SBarry Smith         sum  = b[i];
76817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
76917ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
77017ab2063SBarry Smith       }
77117ab2063SBarry Smith     }
77217ab2063SBarry Smith   }
77317ab2063SBarry Smith   return 0;
77417ab2063SBarry Smith }
77517ab2063SBarry Smith 
776d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
77717ab2063SBarry Smith {
778416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
779416022c9SBarry Smith   *nz      = a->nz;
780416022c9SBarry Smith   *nzalloc = a->maxnz;
781416022c9SBarry Smith   *mem     = (int)A->mem;
78217ab2063SBarry Smith   return 0;
78317ab2063SBarry Smith }
78417ab2063SBarry Smith 
78517ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
78617ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
78717ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
78817ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
78917ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
79017ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
79117ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
79217ab2063SBarry Smith 
79317ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
79417ab2063SBarry Smith {
795416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
796416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
79717ab2063SBarry Smith 
79817ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
79917ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
80017ab2063SBarry Smith   if (diag) {
80117ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
802416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
803416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
804416022c9SBarry Smith         a->ilen[rows[i]] = 1;
805416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
806416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
80717ab2063SBarry Smith       }
80817ab2063SBarry Smith       else {
80917ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
81017ab2063SBarry Smith         CHKERRQ(ierr);
81117ab2063SBarry Smith       }
81217ab2063SBarry Smith     }
81317ab2063SBarry Smith   }
81417ab2063SBarry Smith   else {
81517ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
816416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
817416022c9SBarry Smith       a->ilen[rows[i]] = 0;
81817ab2063SBarry Smith     }
81917ab2063SBarry Smith   }
82017ab2063SBarry Smith   ISRestoreIndices(is,&rows);
82117ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82217ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82317ab2063SBarry Smith   return 0;
82417ab2063SBarry Smith }
82517ab2063SBarry Smith 
826416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
82717ab2063SBarry Smith {
828416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
829416022c9SBarry Smith   *m = a->m; *n = a->n;
83017ab2063SBarry Smith   return 0;
83117ab2063SBarry Smith }
83217ab2063SBarry Smith 
833416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
83417ab2063SBarry Smith {
835416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
836416022c9SBarry Smith   *m = 0; *n = a->m;
83717ab2063SBarry Smith   return 0;
83817ab2063SBarry Smith }
839416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
84017ab2063SBarry Smith {
841416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
842416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
84317ab2063SBarry Smith 
844416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
84517ab2063SBarry Smith 
846416022c9SBarry Smith   if (!a->assembled) {
847416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
848416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
84917ab2063SBarry Smith   }
850416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
851416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
85217ab2063SBarry Smith   if (idx) {
85317ab2063SBarry Smith     if (*nz) {
854416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
8550452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
85617ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
85717ab2063SBarry Smith     }
85817ab2063SBarry Smith     else *idx = 0;
85917ab2063SBarry Smith   }
86017ab2063SBarry Smith   return 0;
86117ab2063SBarry Smith }
86217ab2063SBarry Smith 
863416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
86417ab2063SBarry Smith {
8650452661fSBarry Smith   if (idx) {if (*idx) PetscFree(*idx);}
86617ab2063SBarry Smith   return 0;
86717ab2063SBarry Smith }
86817ab2063SBarry Smith 
869cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
87017ab2063SBarry Smith {
871416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
872416022c9SBarry Smith   Scalar     *v = a->a;
87317ab2063SBarry Smith   double     sum = 0.0;
874416022c9SBarry Smith   int        i, j,shift = a->indexshift;
87517ab2063SBarry Smith 
876416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
87717ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
878416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
87917ab2063SBarry Smith #if defined(PETSC_COMPLEX)
88017ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
88117ab2063SBarry Smith #else
88217ab2063SBarry Smith       sum += (*v)*(*v); v++;
88317ab2063SBarry Smith #endif
88417ab2063SBarry Smith     }
88517ab2063SBarry Smith     *norm = sqrt(sum);
88617ab2063SBarry Smith   }
88717ab2063SBarry Smith   else if (type == NORM_1) {
88817ab2063SBarry Smith     double *tmp;
889416022c9SBarry Smith     int    *jj = a->j;
8900452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
891cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
89217ab2063SBarry Smith     *norm = 0.0;
893416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
894a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
89517ab2063SBarry Smith     }
896416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
89717ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
89817ab2063SBarry Smith     }
8990452661fSBarry Smith     PetscFree(tmp);
90017ab2063SBarry Smith   }
90117ab2063SBarry Smith   else if (type == NORM_INFINITY) {
90217ab2063SBarry Smith     *norm = 0.0;
903416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
904416022c9SBarry Smith       v = a->a + a->i[j] + shift;
90517ab2063SBarry Smith       sum = 0.0;
906416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
907cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
90817ab2063SBarry Smith       }
90917ab2063SBarry Smith       if (sum > *norm) *norm = sum;
91017ab2063SBarry Smith     }
91117ab2063SBarry Smith   }
91217ab2063SBarry Smith   else {
91348d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
91417ab2063SBarry Smith   }
91517ab2063SBarry Smith   return 0;
91617ab2063SBarry Smith }
91717ab2063SBarry Smith 
918416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
91917ab2063SBarry Smith {
920416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
921416022c9SBarry Smith   Mat        C;
922416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
923416022c9SBarry Smith   int        shift = a->indexshift;
924d5d45c9bSBarry Smith   Scalar     *array = a->a;
92517ab2063SBarry Smith 
926416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
9270452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
928cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
92917ab2063SBarry Smith   if (shift) {
93017ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
93117ab2063SBarry Smith   }
93217ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
933416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
9340452661fSBarry Smith   PetscFree(col);
93517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
93617ab2063SBarry Smith     len = ai[i+1]-ai[i];
937416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
93817ab2063SBarry Smith     array += len; aj += len;
93917ab2063SBarry Smith   }
94017ab2063SBarry Smith   if (shift) {
94117ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
94217ab2063SBarry Smith   }
94317ab2063SBarry Smith 
944416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
945416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
94617ab2063SBarry Smith 
947416022c9SBarry Smith   if (B) {
948416022c9SBarry Smith     *B = C;
94917ab2063SBarry Smith   } else {
950416022c9SBarry Smith     /* This isn't really an in-place transpose */
9510452661fSBarry Smith     PetscFree(a->a);
9520452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
9530452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
9540452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
9550452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
9560452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
9570452661fSBarry Smith     PetscFree(a);
958416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
9590452661fSBarry Smith     PetscHeaderDestroy(C);
96017ab2063SBarry Smith   }
96117ab2063SBarry Smith   return 0;
96217ab2063SBarry Smith }
96317ab2063SBarry Smith 
964416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
96517ab2063SBarry Smith {
966416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
96717ab2063SBarry Smith   Scalar     *l,*r,x,*v;
968d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
96917ab2063SBarry Smith 
97048d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
97117ab2063SBarry Smith   if (ll) {
97217ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
973416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
974416022c9SBarry Smith     v = a->a;
97517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
97617ab2063SBarry Smith       x = l[i];
977416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
97817ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
97917ab2063SBarry Smith     }
98017ab2063SBarry Smith   }
98117ab2063SBarry Smith   if (rr) {
98217ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
983416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
984416022c9SBarry Smith     v = a->a; jj = a->j;
98517ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
98617ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
98717ab2063SBarry Smith     }
98817ab2063SBarry Smith   }
98917ab2063SBarry Smith   return 0;
99017ab2063SBarry Smith }
99117ab2063SBarry Smith 
992cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
99317ab2063SBarry Smith {
994db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
99502834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
996a2744918SBarry Smith   register int sum,lensi;
99702834360SBarry Smith   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
998a2744918SBarry Smith   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
999db02288aSLois Curfman McInnes   Scalar       *vwork,*a_new;
1000416022c9SBarry Smith   Mat          C;
100117ab2063SBarry Smith 
1002416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
100317ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
100417ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
100517ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
100617ab2063SBarry Smith 
100702834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
100802834360SBarry Smith     /* special case of contiguous rows */
10090452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
101002834360SBarry Smith     starts = lens + ncols;
101102834360SBarry Smith     /* loop over new rows determining lens and starting points */
101202834360SBarry Smith     for (i=0; i<nrows; i++) {
1013a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1014a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
101502834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1016a2744918SBarry Smith         if (aj[k] >= first) {
101702834360SBarry Smith           starts[i] = k;
101802834360SBarry Smith           break;
101902834360SBarry Smith 	}
102002834360SBarry Smith       }
1021a2744918SBarry Smith       sum = 0;
102202834360SBarry Smith       while (k < kend) {
1023a2744918SBarry Smith         if (aj[k++] >= first+ncols) break;
1024a2744918SBarry Smith         sum++;
102502834360SBarry Smith       }
1026a2744918SBarry Smith       lens[i] = sum;
102702834360SBarry Smith     }
102802834360SBarry Smith     /* create submatrix */
1029cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
103008480c60SBarry Smith       int n_cols,n_rows;
103108480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
103208480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
103308480c60SBarry Smith       MatZeroEntries(*B);
103408480c60SBarry Smith       C = *B;
103508480c60SBarry Smith     }
103608480c60SBarry Smith     else {
103702834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
103808480c60SBarry Smith     }
1039db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1040db02288aSLois Curfman McInnes 
104102834360SBarry Smith     /* loop over rows inserting into submatrix */
1042db02288aSLois Curfman McInnes     a_new    = c->a;
1043db02288aSLois Curfman McInnes     j_new    = c->j;
1044db02288aSLois Curfman McInnes     i_new    = c->i;
1045db02288aSLois Curfman McInnes     i_new[0] = -shift;
104602834360SBarry Smith     for (i=0; i<nrows; i++) {
1047a2744918SBarry Smith       ii    = starts[i];
1048a2744918SBarry Smith       lensi = lens[i];
1049a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1050a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
105102834360SBarry Smith       }
1052a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1053a2744918SBarry Smith       a_new      += lensi;
1054a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1055a2744918SBarry Smith       c->ilen[i]  = lensi;
105602834360SBarry Smith     }
10570452661fSBarry Smith     PetscFree(lens);
105802834360SBarry Smith   }
105902834360SBarry Smith   else {
106002834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10610452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
106202834360SBarry Smith     ssmap = smap + shift;
10637eb43aa7SLois Curfman McInnes     cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork);
106402834360SBarry Smith     lens  = cwork + ncols;
10650452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1066cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
106717ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
106802834360SBarry Smith     /* determine lens of each row */
106902834360SBarry Smith     for (i=0; i<nrows; i++) {
107002834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
107102834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
107202834360SBarry Smith       lens[i] = 0;
107302834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
107402834360SBarry Smith         if (ssmap[a->j[k]]) {
107502834360SBarry Smith           lens[i]++;
107602834360SBarry Smith         }
107702834360SBarry Smith       }
107802834360SBarry Smith     }
107917ab2063SBarry Smith     /* Create and fill new matrix */
1080a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
108108480c60SBarry Smith       int n_cols,n_rows;
108208480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
108308480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
108408480c60SBarry Smith       MatZeroEntries(*B);
108508480c60SBarry Smith       C = *B;
108608480c60SBarry Smith     }
108708480c60SBarry Smith     else {
108802834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
108908480c60SBarry Smith     }
109017ab2063SBarry Smith     for (i=0; i<nrows; i++) {
109117ab2063SBarry Smith       nznew  = 0;
1092416022c9SBarry Smith       kstart = a->i[irow[i]]+shift;
1093416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
109417ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
109502834360SBarry Smith         if (ssmap[a->j[k]]) {
109602834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1097416022c9SBarry Smith           vwork[nznew++] = a->a[k];
109817ab2063SBarry Smith         }
109917ab2063SBarry Smith       }
1100416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
110117ab2063SBarry Smith     }
110202834360SBarry Smith     /* Free work space */
110302834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
11040452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
110502834360SBarry Smith   }
1106416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1107416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
110817ab2063SBarry Smith 
110917ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1110416022c9SBarry Smith   *B = C;
111117ab2063SBarry Smith   return 0;
111217ab2063SBarry Smith }
111317ab2063SBarry Smith 
1114a871dcd8SBarry Smith /*
111563b91edcSBarry Smith      note: This can only work for identity for row and col. It would
111663b91edcSBarry Smith    be good to check this and otherwise generate an error.
1117a871dcd8SBarry Smith */
111863b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1119a871dcd8SBarry Smith {
112063b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
112108480c60SBarry Smith   int        ierr;
112263b91edcSBarry Smith   Mat        outA;
112363b91edcSBarry Smith 
1124a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1125a871dcd8SBarry Smith 
112663b91edcSBarry Smith   outA          = inA;
112763b91edcSBarry Smith   inA->factor   = FACTOR_LU;
112863b91edcSBarry Smith   a->row        = row;
112963b91edcSBarry Smith   a->col        = col;
113063b91edcSBarry Smith 
11310452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
113263b91edcSBarry Smith 
113308480c60SBarry Smith   if (!a->diag) {
113408480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
113563b91edcSBarry Smith   }
113663b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1137a871dcd8SBarry Smith   return 0;
1138a871dcd8SBarry Smith }
1139a871dcd8SBarry Smith 
1140cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1141cddf8d76SBarry Smith                                     Mat **B)
1142cddf8d76SBarry Smith {
1143cddf8d76SBarry Smith   int ierr,i;
1144cddf8d76SBarry Smith 
1145cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
11460452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1147cddf8d76SBarry Smith   }
1148cddf8d76SBarry Smith 
1149cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1150cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1151cddf8d76SBarry Smith   }
1152cddf8d76SBarry Smith   return 0;
1153cddf8d76SBarry Smith }
1154cddf8d76SBarry Smith 
11554dcbc457SBarry Smith static int MatIncreaseOverlap_SeqAIJ(Mat A, int n, IS *is, int ov)
11564dcbc457SBarry Smith {
11574dcbc457SBarry Smith   int i,m,*idx,ierr;
11584dcbc457SBarry Smith 
11594dcbc457SBarry Smith   for ( i=0; i<n; i++ ) {
11604dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr);
11614dcbc457SBarry Smith     ISGetLocalSize(is[i],&m);
11624dcbc457SBarry Smith   }
11634dcbc457SBarry Smith   SETERRQ(1,"MatIncreaseOverlap_SeqAIJ:Not implemented");
11644dcbc457SBarry Smith }
116517ab2063SBarry Smith /* -------------------------------------------------------------------*/
116617ab2063SBarry Smith 
116717ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
116817ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1169416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1170416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
117117ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
117217ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
117317ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
117417ab2063SBarry Smith        MatRelax_SeqAIJ,
117517ab2063SBarry Smith        MatTranspose_SeqAIJ,
117617ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
117717ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
117817ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
117917ab2063SBarry Smith        MatCompress_SeqAIJ,
118017ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
118117ab2063SBarry Smith        MatGetReordering_SeqAIJ,
118217ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
118317ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
118417ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
118517ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1186416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
1187*3d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1188cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
11897eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
11907eb43aa7SLois Curfman McInnes        MatGetValues_SeqAIJ};
119117ab2063SBarry Smith 
119217ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
119317ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
119417ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
119517ab2063SBarry Smith 
119617ab2063SBarry Smith /*@C
119717ab2063SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
119817ab2063SBarry Smith    (the default uniprocessor PETSc format).
119917ab2063SBarry Smith 
120017ab2063SBarry Smith    Input Parameters:
120117ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
120217ab2063SBarry Smith .  m - number of rows
120317ab2063SBarry Smith .  n - number of columns
120417ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
120517ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
120617ab2063SBarry Smith 
120717ab2063SBarry Smith    Output Parameter:
1208416022c9SBarry Smith .  A - the matrix
120917ab2063SBarry Smith 
121017ab2063SBarry Smith    Notes:
121117ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
121217ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
12130002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
12140002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
121517ab2063SBarry Smith 
121617ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1217b4fd4287SBarry Smith    Set nz=0 and nnz=PETSC_NULL for PETSc to control dynamic memory
121817ab2063SBarry Smith    allocation.
121917ab2063SBarry Smith 
12206c7ebb05SLois Curfman McInnes    By default, this format uses inodes (identical nodes) when possible.
12216c7ebb05SLois Curfman McInnes    We search for consecutive rows with the same nonzero structure, thereby
12226c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
12236c7ebb05SLois Curfman McInnes 
12246c7ebb05SLois Curfman McInnes    Options Database Keys:
12256c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
12266c7ebb05SLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)
122717ab2063SBarry Smith 
122817ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
122917ab2063SBarry Smith @*/
1230416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
123117ab2063SBarry Smith {
1232416022c9SBarry Smith   Mat        B;
1233416022c9SBarry Smith   Mat_SeqAIJ *b;
123417ab2063SBarry Smith   int        i,len,ierr;
1235d5d45c9bSBarry Smith 
1236416022c9SBarry Smith   *A      = 0;
12370452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1238416022c9SBarry Smith   PLogObjectCreate(B);
12390452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1240416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1241416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1242416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1243416022c9SBarry Smith   B->factor           = 0;
1244416022c9SBarry Smith   B->lupivotthreshold = 1.0;
1245b4fd4287SBarry Smith   OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1246416022c9SBarry Smith   b->row              = 0;
1247416022c9SBarry Smith   b->col              = 0;
1248416022c9SBarry Smith   b->indexshift       = 0;
1249b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_oneindex")) b->indexshift = -1;
125017ab2063SBarry Smith 
1251416022c9SBarry Smith   b->m       = m;
1252416022c9SBarry Smith   b->n       = n;
12530452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1254b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
125517ab2063SBarry Smith     if (nz <= 0) nz = 1;
1256416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
125717ab2063SBarry Smith     nz = nz*m;
125817ab2063SBarry Smith   }
125917ab2063SBarry Smith   else {
126017ab2063SBarry Smith     nz = 0;
1261416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
126217ab2063SBarry Smith   }
126317ab2063SBarry Smith 
126417ab2063SBarry Smith   /* allocate the matrix space */
1265416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
12660452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1267416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1268cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1269416022c9SBarry Smith   b->i  = b->j + nz;
1270416022c9SBarry Smith   b->singlemalloc = 1;
127117ab2063SBarry Smith 
1272416022c9SBarry Smith   b->i[0] = -b->indexshift;
127317ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1274416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
127517ab2063SBarry Smith   }
127617ab2063SBarry Smith 
1277416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
12780452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1279416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1280416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
128117ab2063SBarry Smith 
1282416022c9SBarry Smith   b->nz               = 0;
1283416022c9SBarry Smith   b->maxnz            = nz;
1284416022c9SBarry Smith   b->sorted           = 0;
1285416022c9SBarry Smith   b->roworiented      = 1;
1286416022c9SBarry Smith   b->nonew            = 0;
1287416022c9SBarry Smith   b->diag             = 0;
1288416022c9SBarry Smith   b->assembled        = 0;
1289416022c9SBarry Smith   b->solve_work       = 0;
129071bd300dSLois Curfman McInnes   b->spptr            = 0;
1291754ec7b1SSatish Balay   b->inode.node_count = 0;
1292754ec7b1SSatish Balay   b->inode.size       = 0;
12936c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
12946c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
129517ab2063SBarry Smith 
1296416022c9SBarry Smith   *A = B;
1297b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_superlu")) {
1298416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
129917ab2063SBarry Smith   }
1300b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_essl")) {
1301416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
130217ab2063SBarry Smith   }
1303b4fd4287SBarry Smith   if (OptionsHasName(PETSC_NULL,"-mat_aij_dxml")) {
1304416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1305416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
130617ab2063SBarry Smith   }
130717ab2063SBarry Smith 
130817ab2063SBarry Smith   return 0;
130917ab2063SBarry Smith }
131017ab2063SBarry Smith 
1311*3d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
131217ab2063SBarry Smith {
1313416022c9SBarry Smith   Mat        C;
1314416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
131508480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
131617ab2063SBarry Smith 
13174043dd9cSLois Curfman McInnes   *B = 0;
1318*3d1612f7SBarry Smith   if (!a->assembled) SETERRQ(1,"MatConvertSameType_SeqAIJ:Cannot copy unassembled matrix");
13190452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1320416022c9SBarry Smith   PLogObjectCreate(C);
13210452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
132241c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1323416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1324416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1325416022c9SBarry Smith   C->factor     = A->factor;
1326416022c9SBarry Smith   c->row        = 0;
1327416022c9SBarry Smith   c->col        = 0;
1328416022c9SBarry Smith   c->indexshift = shift;
132917ab2063SBarry Smith 
1330416022c9SBarry Smith   c->m          = a->m;
1331416022c9SBarry Smith   c->n          = a->n;
133217ab2063SBarry Smith 
13330452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
13340452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
133517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1336416022c9SBarry Smith     c->imax[i] = a->imax[i];
1337416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
133817ab2063SBarry Smith   }
133917ab2063SBarry Smith 
134017ab2063SBarry Smith   /* allocate the matrix space */
1341416022c9SBarry Smith   c->singlemalloc = 1;
1342416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
13430452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1344416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1345416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1346416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
134717ab2063SBarry Smith   if (m > 0) {
1348416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
134908480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1350416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
135117ab2063SBarry Smith     }
135208480c60SBarry Smith   }
135317ab2063SBarry Smith 
1354416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1355416022c9SBarry Smith   c->sorted      = a->sorted;
1356416022c9SBarry Smith   c->roworiented = a->roworiented;
1357416022c9SBarry Smith   c->nonew       = a->nonew;
135817ab2063SBarry Smith 
1359416022c9SBarry Smith   if (a->diag) {
13600452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1361416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
136217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1363416022c9SBarry Smith       c->diag[i] = a->diag[i];
136417ab2063SBarry Smith     }
136517ab2063SBarry Smith   }
1366416022c9SBarry Smith   else c->diag          = 0;
13676c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
13686c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1369754ec7b1SSatish Balay   if (a->inode.size){
1370754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1371754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1372754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1373754ec7b1SSatish Balay   } else {
1374754ec7b1SSatish Balay     c->inode.size       = 0;
1375754ec7b1SSatish Balay     c->inode.node_count = 0;
1376754ec7b1SSatish Balay   }
1377416022c9SBarry Smith   c->assembled          = 1;
1378416022c9SBarry Smith   c->nz                 = a->nz;
1379416022c9SBarry Smith   c->maxnz              = a->maxnz;
1380416022c9SBarry Smith   c->solve_work         = 0;
138176dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1382754ec7b1SSatish Balay 
1383416022c9SBarry Smith   *B = C;
138417ab2063SBarry Smith   return 0;
138517ab2063SBarry Smith }
138617ab2063SBarry Smith 
1387416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
138817ab2063SBarry Smith {
1389416022c9SBarry Smith   Mat_SeqAIJ   *a;
1390416022c9SBarry Smith   Mat          B;
139117699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
139217ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
139317ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
139417ab2063SBarry Smith 
139517699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
139617699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
139717ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1398416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
139948d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
140017ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
140117ab2063SBarry Smith 
140217ab2063SBarry Smith   /* read in row lengths */
14030452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1404416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
140517ab2063SBarry Smith 
140617ab2063SBarry Smith   /* create our matrix */
1407416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1408416022c9SBarry Smith   B = *A;
1409416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1410416022c9SBarry Smith   shift = a->indexshift;
141117ab2063SBarry Smith 
141217ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1413416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
141417ab2063SBarry Smith   if (shift) {
141517ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1416416022c9SBarry Smith       a->j[i] += 1;
141717ab2063SBarry Smith     }
141817ab2063SBarry Smith   }
141917ab2063SBarry Smith 
142017ab2063SBarry Smith   /* read in nonzero values */
1421416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
142217ab2063SBarry Smith 
142317ab2063SBarry Smith   /* set matrix "i" values */
1424416022c9SBarry Smith   a->i[0] = -shift;
142517ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1426416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1427416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
142817ab2063SBarry Smith   }
14290452661fSBarry Smith   PetscFree(rowlengths);
143017ab2063SBarry Smith 
1431416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1432416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
143317ab2063SBarry Smith   return 0;
143417ab2063SBarry Smith }
143517ab2063SBarry Smith 
143617ab2063SBarry Smith 
143717ab2063SBarry Smith 
1438