xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 754ec7b11b9bb3f9ec56a12193efbdd4ea5dcf66)
117ab2063SBarry Smith #ifndef lint
2*754ec7b1SSatish Balay static char vcid[] = "$Id: aij.c,v 1.110 1995/11/06 19:21:06 bsmith Exp balay $";
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;
1917ab2063SBarry Smith 
20416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix");
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);
430452661fSBarry Smith   PetscFree(ia); PetscFree(ja);
4417ab2063SBarry Smith   return 0;
4517ab2063SBarry Smith }
4617ab2063SBarry Smith 
47416022c9SBarry Smith #define CHUNKSIZE   10
4817ab2063SBarry Smith 
4917ab2063SBarry Smith /* This version has row oriented v  */
50416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
5117ab2063SBarry Smith {
52416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
53416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
54416022c9SBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen;
55d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
56416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
5717ab2063SBarry Smith 
5817ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
59416022c9SBarry Smith     row  = im[k];
6017ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
61416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
6217ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
6317ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
64416022c9SBarry Smith     low = 0;
6517ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
66416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
67416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
68416022c9SBarry Smith       col = in[l] - shift; value = *v++;
69416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
70416022c9SBarry Smith       while (high-low > 5) {
71416022c9SBarry Smith         t = (low+high)/2;
72416022c9SBarry Smith         if (rp[t] > col) high = t;
73416022c9SBarry Smith         else             low  = t;
7417ab2063SBarry Smith       }
75416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
7617ab2063SBarry Smith         if (rp[i] > col) break;
7717ab2063SBarry Smith         if (rp[i] == col) {
78416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
7917ab2063SBarry Smith           else                  ap[i] = value;
8017ab2063SBarry Smith           goto noinsert;
8117ab2063SBarry Smith         }
8217ab2063SBarry Smith       }
8317ab2063SBarry Smith       if (nonew) goto noinsert;
8417ab2063SBarry Smith       if (nrow >= rmax) {
8517ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
86416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
8717ab2063SBarry Smith         Scalar *new_a;
8817ab2063SBarry Smith 
8917ab2063SBarry Smith         /* malloc new storage space */
90416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
910452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
9217ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
9317ab2063SBarry Smith         new_i   = new_j + new_nz;
9417ab2063SBarry Smith 
9517ab2063SBarry Smith         /* copy over old data into new slots */
9617ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
97416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
98416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
99416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
100416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
10117ab2063SBarry Smith                                                            len*sizeof(int));
102416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
103416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
10417ab2063SBarry Smith                                                            len*sizeof(Scalar));
10517ab2063SBarry Smith         /* free up old matrix storage */
1060452661fSBarry Smith         PetscFree(a->a);
1070452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
108416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
109416022c9SBarry Smith         a->singlemalloc = 1;
11017ab2063SBarry Smith 
11117ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
112416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
113416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
114416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
11517ab2063SBarry Smith       }
116416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
117416022c9SBarry Smith       /* shift up all the later entries in this row */
118416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
11917ab2063SBarry Smith         rp[ii+1] = rp[ii];
12017ab2063SBarry Smith         ap[ii+1] = ap[ii];
12117ab2063SBarry Smith       }
12217ab2063SBarry Smith       rp[i] = col;
12317ab2063SBarry Smith       ap[i] = value;
12417ab2063SBarry Smith       noinsert:;
125416022c9SBarry Smith       low = i + 1;
12617ab2063SBarry Smith     }
12717ab2063SBarry Smith     ailen[row] = nrow;
12817ab2063SBarry Smith   }
12917ab2063SBarry Smith   return 0;
13017ab2063SBarry Smith }
13117ab2063SBarry Smith 
13217ab2063SBarry Smith #include "draw.h"
13317ab2063SBarry Smith #include "pinclude/pviewer.h"
134416022c9SBarry Smith #include "sysio.h"
13517ab2063SBarry Smith 
136416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
13717ab2063SBarry Smith {
138416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
139416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
14017ab2063SBarry Smith 
141416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1420452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
143416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
144416022c9SBarry Smith   col_lens[1] = a->m;
145416022c9SBarry Smith   col_lens[2] = a->n;
146416022c9SBarry Smith   col_lens[3] = a->nz;
147416022c9SBarry Smith 
148416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
149416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
150416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
15117ab2063SBarry Smith   }
152416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
1530452661fSBarry Smith   PetscFree(col_lens);
154416022c9SBarry Smith 
155416022c9SBarry Smith   /* store column indices (zero start index) */
156416022c9SBarry Smith   if (a->indexshift) {
157416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
15817ab2063SBarry Smith   }
159416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
160416022c9SBarry Smith   if (a->indexshift) {
161416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
16217ab2063SBarry Smith   }
163416022c9SBarry Smith 
164416022c9SBarry Smith   /* store nonzero values */
165416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
16617ab2063SBarry Smith   return 0;
16717ab2063SBarry Smith }
168416022c9SBarry Smith 
169416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
170416022c9SBarry Smith {
171416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
172416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
17317ab2063SBarry Smith   FILE        *fd;
17417ab2063SBarry Smith   char        *outputname;
17517ab2063SBarry Smith 
176416022c9SBarry Smith   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
177416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
178416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
17917ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
18008480c60SBarry Smith     /* no need to print additional information */ ;
18117ab2063SBarry Smith   }
18217ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
18317ab2063SBarry Smith     int nz, nzalloc, mem;
184416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
185416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
18617ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
18717ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
18817ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
18917ab2063SBarry Smith 
19017ab2063SBarry Smith     for (i=0; i<m; i++) {
191416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
19217ab2063SBarry Smith #if defined(PETSC_COMPLEX)
193416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
194416022c9SBarry Smith                    imag(a->a[j]));
19517ab2063SBarry Smith #else
196416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
19717ab2063SBarry Smith #endif
19817ab2063SBarry Smith       }
19917ab2063SBarry Smith     }
20017ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
20117ab2063SBarry Smith   }
20217ab2063SBarry Smith   else {
20317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
20417ab2063SBarry Smith       fprintf(fd,"row %d:",i);
205416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
20617ab2063SBarry Smith #if defined(PETSC_COMPLEX)
207416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
208416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
20917ab2063SBarry Smith         }
21017ab2063SBarry Smith         else {
211416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
21217ab2063SBarry Smith         }
21317ab2063SBarry Smith #else
214416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
21517ab2063SBarry Smith #endif
21617ab2063SBarry Smith       }
21717ab2063SBarry Smith       fprintf(fd,"\n");
21817ab2063SBarry Smith     }
21917ab2063SBarry Smith   }
22017ab2063SBarry Smith   fflush(fd);
221416022c9SBarry Smith   return 0;
222416022c9SBarry Smith }
223416022c9SBarry Smith 
224416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
225416022c9SBarry Smith {
226416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
227cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
228cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
229416022c9SBarry Smith   DrawCtx     draw = (DrawCtx) viewer;
230cddf8d76SBarry Smith   DrawButton  button;
231cddf8d76SBarry Smith 
232416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
233416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
234416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
235416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
236cddf8d76SBarry Smith   color = DRAW_BLUE;
237416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
238cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
239416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
240cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
241cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
242cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
243cddf8d76SBarry Smith #else
244cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
245cddf8d76SBarry Smith #endif
246cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
247cddf8d76SBarry Smith     }
248cddf8d76SBarry Smith   }
249cddf8d76SBarry Smith   color = DRAW_CYAN;
250cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
251cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
252cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
253cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
254cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
255cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
256cddf8d76SBarry Smith     }
257cddf8d76SBarry Smith   }
258cddf8d76SBarry Smith   color = DRAW_RED;
259cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
260cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
261cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
262cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
263cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
264cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
265cddf8d76SBarry Smith #else
266cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
267cddf8d76SBarry Smith #endif
268cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
269416022c9SBarry Smith     }
270416022c9SBarry Smith   }
271416022c9SBarry Smith   DrawFlush(draw);
272cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
273cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
274cddf8d76SBarry Smith 
275cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
276cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
277cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
278cddf8d76SBarry Smith     DrawClear(draw);
279cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
280cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
281cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
282cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
283cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
284cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
285cddf8d76SBarry Smith     w *= scale; h *= scale;
286cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
287cddf8d76SBarry Smith     color = DRAW_BLUE;
288cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
289cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
290cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
291cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
292cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
293cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
294cddf8d76SBarry Smith #else
295cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
296cddf8d76SBarry Smith #endif
297cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
298cddf8d76SBarry Smith       }
299cddf8d76SBarry Smith     }
300cddf8d76SBarry Smith     color = DRAW_CYAN;
301cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
302cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
303cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
304cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
305cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
306cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
307cddf8d76SBarry Smith       }
308cddf8d76SBarry Smith     }
309cddf8d76SBarry Smith     color = DRAW_RED;
310cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
311cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
312cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
313cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
314cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
315cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
316cddf8d76SBarry Smith #else
317cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
318cddf8d76SBarry Smith #endif
319cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
320cddf8d76SBarry Smith       }
321cddf8d76SBarry Smith     }
322cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
323cddf8d76SBarry Smith   }
324416022c9SBarry Smith   return 0;
325416022c9SBarry Smith }
326416022c9SBarry Smith 
327416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
328416022c9SBarry Smith {
329416022c9SBarry Smith   Mat         A = (Mat) obj;
330416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
331416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
332416022c9SBarry Smith 
333416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
334416022c9SBarry Smith   if (!viewer) {
335416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
336416022c9SBarry Smith   }
337416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
338416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
339416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
340416022c9SBarry Smith     }
341416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
342416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
343416022c9SBarry Smith     }
344416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
345416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
346416022c9SBarry Smith     }
347416022c9SBarry Smith   }
348416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
349416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
350416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
35117ab2063SBarry Smith   }
35217ab2063SBarry Smith   return 0;
35317ab2063SBarry Smith }
35417ab2063SBarry Smith 
355416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
35617ab2063SBarry Smith {
357416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
358416022c9SBarry Smith   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
359416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
360416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
36117ab2063SBarry Smith 
36217ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
36317ab2063SBarry Smith 
36417ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
365416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
36617ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
36717ab2063SBarry Smith     if (fshift) {
368416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
36917ab2063SBarry Smith       N = ailen[i];
37017ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
37117ab2063SBarry Smith         ip[j-fshift] = ip[j];
37217ab2063SBarry Smith         ap[j-fshift] = ap[j];
37317ab2063SBarry Smith       }
37417ab2063SBarry Smith     }
37517ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
37617ab2063SBarry Smith   }
37717ab2063SBarry Smith   if (m) {
37817ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
37917ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
38017ab2063SBarry Smith   }
38117ab2063SBarry Smith   /* reset ilen and imax for each row */
38217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
38317ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
38417ab2063SBarry Smith   }
385416022c9SBarry Smith   a->nz = ai[m] + shift;
38617ab2063SBarry Smith 
38717ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
388416022c9SBarry Smith   if (fshift && a->diag) {
3890452661fSBarry Smith     PetscFree(a->diag);
390416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
391416022c9SBarry Smith     a->diag = 0;
39217ab2063SBarry Smith   }
393416022c9SBarry Smith   a->assembled = 1;
39417ab2063SBarry Smith   return 0;
39517ab2063SBarry Smith }
39617ab2063SBarry Smith 
397416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
39817ab2063SBarry Smith {
399416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
400cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
40117ab2063SBarry Smith   return 0;
40217ab2063SBarry Smith }
403416022c9SBarry Smith 
40417ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
40517ab2063SBarry Smith {
406416022c9SBarry Smith   Mat        A  = (Mat) obj;
407416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
408d5d45c9bSBarry Smith 
40917ab2063SBarry Smith #if defined(PETSC_LOG)
410416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
41117ab2063SBarry Smith #endif
4120452661fSBarry Smith   PetscFree(a->a);
4130452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4140452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
4150452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4160452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
4170452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
4180452661fSBarry Smith   PetscFree(a);
419416022c9SBarry Smith   PLogObjectDestroy(A);
4200452661fSBarry Smith   PetscHeaderDestroy(A);
42117ab2063SBarry Smith   return 0;
42217ab2063SBarry Smith }
42317ab2063SBarry Smith 
424416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
42517ab2063SBarry Smith {
42617ab2063SBarry Smith   return 0;
42717ab2063SBarry Smith }
42817ab2063SBarry Smith 
429416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
43017ab2063SBarry Smith {
431416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
432416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
433416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
434416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
435416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
436e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
437e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
438e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
439e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
440df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
441df719601SLois Curfman McInnes   else if (op == COLUMN_ORIENTED)
4424d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");}
443df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
4444d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
445e2f28af5SBarry Smith   else
4464d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
44717ab2063SBarry Smith   return 0;
44817ab2063SBarry Smith }
44917ab2063SBarry Smith 
450416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
45117ab2063SBarry Smith {
452416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
453416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
45417ab2063SBarry Smith   Scalar     *x, zero = 0.0;
45517ab2063SBarry Smith 
456416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
45717ab2063SBarry Smith   VecSet(&zero,v);
45817ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
459416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
460416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
461416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
462416022c9SBarry Smith       if (a->j[j]+shift == i) {
463416022c9SBarry Smith         x[i] = a->a[j];
46417ab2063SBarry Smith         break;
46517ab2063SBarry Smith       }
46617ab2063SBarry Smith     }
46717ab2063SBarry Smith   }
46817ab2063SBarry Smith   return 0;
46917ab2063SBarry Smith }
47017ab2063SBarry Smith 
47117ab2063SBarry Smith /* -------------------------------------------------------*/
47217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
47317ab2063SBarry Smith /* -------------------------------------------------------*/
474416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
47517ab2063SBarry Smith {
476416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
47717ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
478416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
47917ab2063SBarry Smith 
480416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
48117ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
482cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
48317ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
48417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
485416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
486416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
487416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
48817ab2063SBarry Smith     alpha = x[i];
48917ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
49017ab2063SBarry Smith   }
491416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
49217ab2063SBarry Smith   return 0;
49317ab2063SBarry Smith }
494d5d45c9bSBarry Smith 
495416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
49617ab2063SBarry Smith {
497416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
49817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
499416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
50017ab2063SBarry Smith 
501416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
50217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
50317ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
50417ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
50517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
506416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
507416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
508416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
50917ab2063SBarry Smith     alpha = x[i];
51017ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
51117ab2063SBarry Smith   }
51217ab2063SBarry Smith   return 0;
51317ab2063SBarry Smith }
51417ab2063SBarry Smith 
515416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
51617ab2063SBarry Smith {
517416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
51817ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
519416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
52017ab2063SBarry Smith 
521416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
52217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
52317ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
524416022c9SBarry Smith   idx  = a->j;
525416022c9SBarry Smith   v    = a->a;
526416022c9SBarry Smith   ii   = a->i;
52717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
528416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
52917ab2063SBarry Smith     sum  = 0.0;
53017ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
53117ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
532416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
53317ab2063SBarry Smith     y[i] = sum;
53417ab2063SBarry Smith   }
535416022c9SBarry Smith   PLogFlops(2*a->nz - m);
53617ab2063SBarry Smith   return 0;
53717ab2063SBarry Smith }
53817ab2063SBarry Smith 
539416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
54017ab2063SBarry Smith {
541416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
54217ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
543cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
54417ab2063SBarry Smith 
54548d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
54617ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
54717ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
548cddf8d76SBarry Smith   idx  = a->j;
549cddf8d76SBarry Smith   v    = a->a;
550cddf8d76SBarry Smith   ii   = a->i;
55117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
552cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
55317ab2063SBarry Smith     sum  = y[i];
554cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
555cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
55617ab2063SBarry Smith     z[i] = sum;
55717ab2063SBarry Smith   }
558416022c9SBarry Smith   PLogFlops(2*a->nz);
55917ab2063SBarry Smith   return 0;
56017ab2063SBarry Smith }
56117ab2063SBarry Smith 
56217ab2063SBarry Smith /*
56317ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
56417ab2063SBarry Smith */
56517ab2063SBarry Smith 
566416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
56717ab2063SBarry Smith {
568416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
569416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
57017ab2063SBarry Smith 
571416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
5720452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
573416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
574416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
575416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
576416022c9SBarry Smith       if (a->j[j]+shift == i) {
57717ab2063SBarry Smith         diag[i] = j - shift;
57817ab2063SBarry Smith         break;
57917ab2063SBarry Smith       }
58017ab2063SBarry Smith     }
58117ab2063SBarry Smith   }
582416022c9SBarry Smith   a->diag = diag;
58317ab2063SBarry Smith   return 0;
58417ab2063SBarry Smith }
58517ab2063SBarry Smith 
586416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
58717ab2063SBarry Smith                            double fshift,int its,Vec xx)
58817ab2063SBarry Smith {
589416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
590416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
591d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
59217ab2063SBarry Smith 
59317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
594416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
595416022c9SBarry Smith   diag = a->diag;
596416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
59717ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
59817ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
59917ab2063SBarry Smith     bs = b + shift;
60017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
601416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
602416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
603416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
604416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
60517ab2063SBarry Smith         sum  = b[i]*d/omega;
60617ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
60717ab2063SBarry Smith         x[i] = sum;
60817ab2063SBarry Smith     }
60917ab2063SBarry Smith     return 0;
61017ab2063SBarry Smith   }
61117ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
61217ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
61317ab2063SBarry Smith   }
614416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
61517ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
61617ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
61717ab2063SBarry Smith 
61817ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
61917ab2063SBarry Smith 
62017ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
62117ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
62217ab2063SBarry Smith     is the relaxation factor.
62317ab2063SBarry Smith     */
6240452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
62517ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
62617ab2063SBarry Smith 
62717ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
62817ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
629416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
630416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
631416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
632416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
63317ab2063SBarry Smith       sum  = b[i];
63417ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
63517ab2063SBarry Smith       x[i] = omega*(sum/d);
63617ab2063SBarry Smith     }
63717ab2063SBarry Smith 
63817ab2063SBarry Smith     /*  t = b - (2*E - D)x */
639416022c9SBarry Smith     v = a->a;
64017ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
64117ab2063SBarry Smith 
64217ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
643416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
644416022c9SBarry Smith     diag = a->diag;
64517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
646416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
647416022c9SBarry Smith       n    = diag[i] - a->i[i];
648416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
649416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
65017ab2063SBarry Smith       sum  = t[i];
65117ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
65217ab2063SBarry Smith       t[i] = omega*(sum/d);
65317ab2063SBarry Smith     }
65417ab2063SBarry Smith 
65517ab2063SBarry Smith     /*  x = x + t */
65617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
6570452661fSBarry Smith     PetscFree(t);
65817ab2063SBarry Smith     return 0;
65917ab2063SBarry Smith   }
66017ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
66117ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
66217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
663416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
664416022c9SBarry Smith         n    = diag[i] - a->i[i];
665416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
666416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
66717ab2063SBarry Smith         sum  = b[i];
66817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
66917ab2063SBarry Smith         x[i] = omega*(sum/d);
67017ab2063SBarry Smith       }
67117ab2063SBarry Smith       xb = x;
67217ab2063SBarry Smith     }
67317ab2063SBarry Smith     else xb = b;
67417ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
67517ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
67617ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
677416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
67817ab2063SBarry Smith       }
67917ab2063SBarry Smith     }
68017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
68117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
682416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
683416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
684416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
685416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
68617ab2063SBarry Smith         sum  = xb[i];
68717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
68817ab2063SBarry Smith         x[i] = omega*(sum/d);
68917ab2063SBarry Smith       }
69017ab2063SBarry Smith     }
69117ab2063SBarry Smith     its--;
69217ab2063SBarry Smith   }
69317ab2063SBarry Smith   while (its--) {
69417ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
69517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
696416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
697416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
698416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
699416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
70017ab2063SBarry Smith         sum  = b[i];
70117ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
70217ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
70317ab2063SBarry Smith       }
70417ab2063SBarry Smith     }
70517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
70617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
707416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
708416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
709416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
710416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
71117ab2063SBarry Smith         sum  = b[i];
71217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
71317ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
71417ab2063SBarry Smith       }
71517ab2063SBarry Smith     }
71617ab2063SBarry Smith   }
71717ab2063SBarry Smith   return 0;
71817ab2063SBarry Smith }
71917ab2063SBarry Smith 
720d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
72117ab2063SBarry Smith {
722416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
723416022c9SBarry Smith   *nz      = a->nz;
724416022c9SBarry Smith   *nzalloc = a->maxnz;
725416022c9SBarry Smith   *mem     = (int)A->mem;
72617ab2063SBarry Smith   return 0;
72717ab2063SBarry Smith }
72817ab2063SBarry Smith 
72917ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
73017ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
73117ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
73217ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
73317ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
73417ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
73517ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
73617ab2063SBarry Smith 
73717ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
73817ab2063SBarry Smith {
739416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
740416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
74117ab2063SBarry Smith 
74217ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
74317ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
74417ab2063SBarry Smith   if (diag) {
74517ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
746416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
747416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
748416022c9SBarry Smith         a->ilen[rows[i]] = 1;
749416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
750416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
75117ab2063SBarry Smith       }
75217ab2063SBarry Smith       else {
75317ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
75417ab2063SBarry Smith         CHKERRQ(ierr);
75517ab2063SBarry Smith       }
75617ab2063SBarry Smith     }
75717ab2063SBarry Smith   }
75817ab2063SBarry Smith   else {
75917ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
760416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
761416022c9SBarry Smith       a->ilen[rows[i]] = 0;
76217ab2063SBarry Smith     }
76317ab2063SBarry Smith   }
76417ab2063SBarry Smith   ISRestoreIndices(is,&rows);
76517ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
76617ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
76717ab2063SBarry Smith   return 0;
76817ab2063SBarry Smith }
76917ab2063SBarry Smith 
770416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
77117ab2063SBarry Smith {
772416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
773416022c9SBarry Smith   *m = a->m; *n = a->n;
77417ab2063SBarry Smith   return 0;
77517ab2063SBarry Smith }
77617ab2063SBarry Smith 
777416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
77817ab2063SBarry Smith {
779416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
780416022c9SBarry Smith   *m = 0; *n = a->m;
78117ab2063SBarry Smith   return 0;
78217ab2063SBarry Smith }
783416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
78417ab2063SBarry Smith {
785416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
786416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
78717ab2063SBarry Smith 
788416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
78917ab2063SBarry Smith 
790416022c9SBarry Smith   if (!a->assembled) {
791416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
792416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
79317ab2063SBarry Smith   }
794416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
795416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
79617ab2063SBarry Smith   if (idx) {
79717ab2063SBarry Smith     if (*nz) {
798416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
7990452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
80017ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
80117ab2063SBarry Smith     }
80217ab2063SBarry Smith     else *idx = 0;
80317ab2063SBarry Smith   }
80417ab2063SBarry Smith   return 0;
80517ab2063SBarry Smith }
80617ab2063SBarry Smith 
807416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
80817ab2063SBarry Smith {
8090452661fSBarry Smith   if (idx) {if (*idx) PetscFree(*idx);}
81017ab2063SBarry Smith   return 0;
81117ab2063SBarry Smith }
81217ab2063SBarry Smith 
813cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
81417ab2063SBarry Smith {
815416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
816416022c9SBarry Smith   Scalar     *v = a->a;
81717ab2063SBarry Smith   double     sum = 0.0;
818416022c9SBarry Smith   int        i, j,shift = a->indexshift;
81917ab2063SBarry Smith 
820416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
82117ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
822416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
82317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
82417ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
82517ab2063SBarry Smith #else
82617ab2063SBarry Smith       sum += (*v)*(*v); v++;
82717ab2063SBarry Smith #endif
82817ab2063SBarry Smith     }
82917ab2063SBarry Smith     *norm = sqrt(sum);
83017ab2063SBarry Smith   }
83117ab2063SBarry Smith   else if (type == NORM_1) {
83217ab2063SBarry Smith     double *tmp;
833416022c9SBarry Smith     int    *jj = a->j;
8340452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
835cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
83617ab2063SBarry Smith     *norm = 0.0;
837416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
838a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
83917ab2063SBarry Smith     }
840416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
84117ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
84217ab2063SBarry Smith     }
8430452661fSBarry Smith     PetscFree(tmp);
84417ab2063SBarry Smith   }
84517ab2063SBarry Smith   else if (type == NORM_INFINITY) {
84617ab2063SBarry Smith     *norm = 0.0;
847416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
848416022c9SBarry Smith       v = a->a + a->i[j] + shift;
84917ab2063SBarry Smith       sum = 0.0;
850416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
851cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
85217ab2063SBarry Smith       }
85317ab2063SBarry Smith       if (sum > *norm) *norm = sum;
85417ab2063SBarry Smith     }
85517ab2063SBarry Smith   }
85617ab2063SBarry Smith   else {
85748d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
85817ab2063SBarry Smith   }
85917ab2063SBarry Smith   return 0;
86017ab2063SBarry Smith }
86117ab2063SBarry Smith 
862416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
86317ab2063SBarry Smith {
864416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
865416022c9SBarry Smith   Mat        C;
866416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
867416022c9SBarry Smith   int        shift = a->indexshift;
868d5d45c9bSBarry Smith   Scalar     *array = a->a;
86917ab2063SBarry Smith 
870416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
8710452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
872cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
87317ab2063SBarry Smith   if (shift) {
87417ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
87517ab2063SBarry Smith   }
87617ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
877416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
8780452661fSBarry Smith   PetscFree(col);
87917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
88017ab2063SBarry Smith     len = ai[i+1]-ai[i];
881416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
88217ab2063SBarry Smith     array += len; aj += len;
88317ab2063SBarry Smith   }
88417ab2063SBarry Smith   if (shift) {
88517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
88617ab2063SBarry Smith   }
88717ab2063SBarry Smith 
888416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
889416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
89017ab2063SBarry Smith 
891416022c9SBarry Smith   if (B) {
892416022c9SBarry Smith     *B = C;
89317ab2063SBarry Smith   } else {
894416022c9SBarry Smith     /* This isn't really an in-place transpose */
8950452661fSBarry Smith     PetscFree(a->a);
8960452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
8970452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
8980452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
8990452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
9000452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
9010452661fSBarry Smith     PetscFree(a);
902416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
9030452661fSBarry Smith     PetscHeaderDestroy(C);
90417ab2063SBarry Smith   }
90517ab2063SBarry Smith   return 0;
90617ab2063SBarry Smith }
90717ab2063SBarry Smith 
908416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
90917ab2063SBarry Smith {
910416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
91117ab2063SBarry Smith   Scalar     *l,*r,x,*v;
912d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
91317ab2063SBarry Smith 
91448d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
91517ab2063SBarry Smith   if (ll) {
91617ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
917416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
918416022c9SBarry Smith     v = a->a;
91917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
92017ab2063SBarry Smith       x = l[i];
921416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
92217ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
92317ab2063SBarry Smith     }
92417ab2063SBarry Smith   }
92517ab2063SBarry Smith   if (rr) {
92617ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
927416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
928416022c9SBarry Smith     v = a->a; jj = a->j;
92917ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
93017ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
93117ab2063SBarry Smith     }
93217ab2063SBarry Smith   }
93317ab2063SBarry Smith   return 0;
93417ab2063SBarry Smith }
93517ab2063SBarry Smith 
936cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
93717ab2063SBarry Smith {
938db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
93902834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
940a2744918SBarry Smith   register int sum,lensi;
94102834360SBarry Smith   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
942a2744918SBarry Smith   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
943db02288aSLois Curfman McInnes   Scalar       *vwork,*a_new;
944416022c9SBarry Smith   Mat          C;
94517ab2063SBarry Smith 
946416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
94717ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
94817ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
94917ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
95017ab2063SBarry Smith 
95102834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
95202834360SBarry Smith     /* special case of contiguous rows */
9530452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
95402834360SBarry Smith     starts = lens + ncols;
95502834360SBarry Smith     /* loop over new rows determining lens and starting points */
95602834360SBarry Smith     for (i=0; i<nrows; i++) {
957a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
958a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
95902834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
960a2744918SBarry Smith         if (aj[k] >= first) {
96102834360SBarry Smith           starts[i] = k;
96202834360SBarry Smith           break;
96302834360SBarry Smith 	}
96402834360SBarry Smith       }
965a2744918SBarry Smith       sum = 0;
96602834360SBarry Smith       while (k < kend) {
967a2744918SBarry Smith         if (aj[k++] >= first+ncols) break;
968a2744918SBarry Smith         sum++;
96902834360SBarry Smith       }
970a2744918SBarry Smith       lens[i] = sum;
97102834360SBarry Smith     }
97202834360SBarry Smith     /* create submatrix */
973cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
97408480c60SBarry Smith       int n_cols,n_rows;
97508480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
97608480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
97708480c60SBarry Smith       MatZeroEntries(*B);
97808480c60SBarry Smith       C = *B;
97908480c60SBarry Smith     }
98008480c60SBarry Smith     else {
98102834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
98208480c60SBarry Smith     }
983db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
984db02288aSLois Curfman McInnes 
98502834360SBarry Smith     /* loop over rows inserting into submatrix */
986db02288aSLois Curfman McInnes     a_new    = c->a;
987db02288aSLois Curfman McInnes     j_new    = c->j;
988db02288aSLois Curfman McInnes     i_new    = c->i;
989db02288aSLois Curfman McInnes     i_new[0] = -shift;
99002834360SBarry Smith     for (i=0; i<nrows; i++) {
991a2744918SBarry Smith       ii    = starts[i];
992a2744918SBarry Smith       lensi = lens[i];
993a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
994a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
99502834360SBarry Smith       }
996a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
997a2744918SBarry Smith       a_new      += lensi;
998a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
999a2744918SBarry Smith       c->ilen[i]  = lensi;
100002834360SBarry Smith     }
10010452661fSBarry Smith     PetscFree(lens);
100202834360SBarry Smith   }
100302834360SBarry Smith   else {
100402834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10050452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
100602834360SBarry Smith     ssmap = smap + shift;
10070452661fSBarry Smith     cwork = (int *) PetscMalloc((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork);
100802834360SBarry Smith     lens  = cwork + ncols;
10090452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1010cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
101117ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
101202834360SBarry Smith     /* determine lens of each row */
101302834360SBarry Smith     for (i=0; i<nrows; i++) {
101402834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
101502834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
101602834360SBarry Smith       lens[i] = 0;
101702834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
101802834360SBarry Smith         if (ssmap[a->j[k]]) {
101902834360SBarry Smith           lens[i]++;
102002834360SBarry Smith         }
102102834360SBarry Smith       }
102202834360SBarry Smith     }
102317ab2063SBarry Smith     /* Create and fill new matrix */
1024a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
102508480c60SBarry Smith       int n_cols,n_rows;
102608480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
102708480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
102808480c60SBarry Smith       MatZeroEntries(*B);
102908480c60SBarry Smith       C = *B;
103008480c60SBarry Smith     }
103108480c60SBarry Smith     else {
103202834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
103308480c60SBarry Smith     }
103417ab2063SBarry Smith     for (i=0; i<nrows; i++) {
103517ab2063SBarry Smith       nznew  = 0;
1036416022c9SBarry Smith       kstart = a->i[irow[i]]+shift;
1037416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
103817ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
103902834360SBarry Smith         if (ssmap[a->j[k]]) {
104002834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1041416022c9SBarry Smith           vwork[nznew++] = a->a[k];
104217ab2063SBarry Smith         }
104317ab2063SBarry Smith       }
1044416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
104517ab2063SBarry Smith     }
104602834360SBarry Smith     /* Free work space */
104702834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10480452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
104902834360SBarry Smith   }
1050416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1051416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
105217ab2063SBarry Smith 
105317ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1054416022c9SBarry Smith   *B = C;
105517ab2063SBarry Smith   return 0;
105617ab2063SBarry Smith }
105717ab2063SBarry Smith 
1058a871dcd8SBarry Smith /*
105963b91edcSBarry Smith      note: This can only work for identity for row and col. It would
106063b91edcSBarry Smith    be good to check this and otherwise generate an error.
1061a871dcd8SBarry Smith */
106263b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1063a871dcd8SBarry Smith {
106463b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
106508480c60SBarry Smith   int        ierr;
106663b91edcSBarry Smith   Mat        outA;
106763b91edcSBarry Smith 
1068a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1069a871dcd8SBarry Smith 
107063b91edcSBarry Smith   outA          = inA;
107163b91edcSBarry Smith   inA->factor   = FACTOR_LU;
107263b91edcSBarry Smith   a->row        = row;
107363b91edcSBarry Smith   a->col        = col;
107463b91edcSBarry Smith 
10750452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
107663b91edcSBarry Smith 
107708480c60SBarry Smith   if (!a->diag) {
107808480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
107963b91edcSBarry Smith   }
108063b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1081a871dcd8SBarry Smith   return 0;
1082a871dcd8SBarry Smith }
1083a871dcd8SBarry Smith 
1084cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1085cddf8d76SBarry Smith                                     Mat **B)
1086cddf8d76SBarry Smith {
1087cddf8d76SBarry Smith   int ierr,i;
1088cddf8d76SBarry Smith 
1089cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
10900452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1091cddf8d76SBarry Smith   }
1092cddf8d76SBarry Smith 
1093cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1094cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1095cddf8d76SBarry Smith   }
1096cddf8d76SBarry Smith   return 0;
1097cddf8d76SBarry Smith }
1098cddf8d76SBarry Smith 
109917ab2063SBarry Smith /* -------------------------------------------------------------------*/
110017ab2063SBarry Smith 
110117ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
110217ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1103416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1104416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
110517ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
110617ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
110717ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
110817ab2063SBarry Smith        MatRelax_SeqAIJ,
110917ab2063SBarry Smith        MatTranspose_SeqAIJ,
111017ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
111117ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
111217ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
111317ab2063SBarry Smith        MatCompress_SeqAIJ,
111417ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
111517ab2063SBarry Smith        MatGetReordering_SeqAIJ,
111617ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
111717ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
111817ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
111917ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1120416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
1121a871dcd8SBarry Smith        MatCopyPrivate_SeqAIJ,0,0,
1122cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
1123cddf8d76SBarry Smith        MatGetSubMatrices_SeqAIJ};
112417ab2063SBarry Smith 
112517ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
112617ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
112717ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
112817ab2063SBarry Smith 
112917ab2063SBarry Smith /*@C
113017ab2063SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
113117ab2063SBarry Smith    (the default uniprocessor PETSc format).
113217ab2063SBarry Smith 
113317ab2063SBarry Smith    Input Parameters:
113417ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
113517ab2063SBarry Smith .  m - number of rows
113617ab2063SBarry Smith .  n - number of columns
113717ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
113817ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
113917ab2063SBarry Smith 
114017ab2063SBarry Smith    Output Parameter:
1141416022c9SBarry Smith .  A - the matrix
114217ab2063SBarry Smith 
114317ab2063SBarry Smith    Notes:
114417ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
114517ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
11460002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
11470002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
114817ab2063SBarry Smith 
114917ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
115017ab2063SBarry Smith    Set both nz and nnz to zero for PETSc to control dynamic memory
115117ab2063SBarry Smith    allocation.
115217ab2063SBarry Smith 
115317ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse
115417ab2063SBarry Smith 
115517ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
115617ab2063SBarry Smith @*/
1157416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
115817ab2063SBarry Smith {
1159416022c9SBarry Smith   Mat        B;
1160416022c9SBarry Smith   Mat_SeqAIJ *b;
116117ab2063SBarry Smith   int        i,len,ierr;
1162d5d45c9bSBarry Smith 
1163416022c9SBarry Smith   *A      = 0;
11640452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1165416022c9SBarry Smith   PLogObjectCreate(B);
11660452661fSBarry Smith   B->data               = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1167416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1168416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1169416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1170416022c9SBarry Smith   B->factor           = 0;
1171416022c9SBarry Smith   B->lupivotthreshold = 1.0;
1172416022c9SBarry Smith   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1173416022c9SBarry Smith   b->row              = 0;
1174416022c9SBarry Smith   b->col              = 0;
1175416022c9SBarry Smith   b->indexshift       = 0;
1176416022c9SBarry Smith   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
117717ab2063SBarry Smith 
1178416022c9SBarry Smith   b->m       = m;
1179416022c9SBarry Smith   b->n       = n;
11800452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
118117ab2063SBarry Smith   if (!nnz) {
118217ab2063SBarry Smith     if (nz <= 0) nz = 1;
1183416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
118417ab2063SBarry Smith     nz = nz*m;
118517ab2063SBarry Smith   }
118617ab2063SBarry Smith   else {
118717ab2063SBarry Smith     nz = 0;
1188416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
118917ab2063SBarry Smith   }
119017ab2063SBarry Smith 
119117ab2063SBarry Smith   /* allocate the matrix space */
1192416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
11930452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1194416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1195cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1196416022c9SBarry Smith   b->i  = b->j + nz;
1197416022c9SBarry Smith   b->singlemalloc = 1;
119817ab2063SBarry Smith 
1199416022c9SBarry Smith   b->i[0] = -b->indexshift;
120017ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1201416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
120217ab2063SBarry Smith   }
120317ab2063SBarry Smith 
1204416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
12050452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1206416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1207416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
120817ab2063SBarry Smith 
1209416022c9SBarry Smith   b->nz               = 0;
1210416022c9SBarry Smith   b->maxnz            = nz;
1211416022c9SBarry Smith   b->sorted           = 0;
1212416022c9SBarry Smith   b->roworiented      = 1;
1213416022c9SBarry Smith   b->nonew            = 0;
1214416022c9SBarry Smith   b->diag             = 0;
1215416022c9SBarry Smith   b->assembled        = 0;
1216416022c9SBarry Smith   b->solve_work       = 0;
121771bd300dSLois Curfman McInnes   b->spptr            = 0;
1218*754ec7b1SSatish Balay   b->inode.node_count = 0;
1219*754ec7b1SSatish Balay   b->inode.size       = 0;
122017ab2063SBarry Smith 
1221416022c9SBarry Smith   *A = B;
122217ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_superlu")) {
1223416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
122417ab2063SBarry Smith   }
122517ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_essl")) {
1226416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
122717ab2063SBarry Smith   }
122817ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_dxml")) {
1229416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1230416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
123117ab2063SBarry Smith   }
123217ab2063SBarry Smith 
123317ab2063SBarry Smith   return 0;
123417ab2063SBarry Smith }
123517ab2063SBarry Smith 
123608480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues)
123717ab2063SBarry Smith {
1238416022c9SBarry Smith   Mat        C;
1239416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
124008480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
124117ab2063SBarry Smith 
12424043dd9cSLois Curfman McInnes   *B = 0;
1243416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
12440452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1245416022c9SBarry Smith   PLogObjectCreate(C);
12460452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1247416022c9SBarry Smith   PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps));
1248416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1249416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1250416022c9SBarry Smith   C->factor     = A->factor;
1251416022c9SBarry Smith   c->row        = 0;
1252416022c9SBarry Smith   c->col        = 0;
1253416022c9SBarry Smith   c->indexshift = shift;
125417ab2063SBarry Smith 
1255416022c9SBarry Smith   c->m          = a->m;
1256416022c9SBarry Smith   c->n          = a->n;
125717ab2063SBarry Smith 
12580452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
12590452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
126017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1261416022c9SBarry Smith     c->imax[i] = a->imax[i];
1262416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
126317ab2063SBarry Smith   }
126417ab2063SBarry Smith 
126517ab2063SBarry Smith   /* allocate the matrix space */
1266416022c9SBarry Smith   c->singlemalloc = 1;
1267416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
12680452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1269416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1270416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1271416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
127217ab2063SBarry Smith   if (m > 0) {
1273416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
127408480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1275416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
127617ab2063SBarry Smith     }
127708480c60SBarry Smith   }
127817ab2063SBarry Smith 
1279416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1280416022c9SBarry Smith   c->sorted      = a->sorted;
1281416022c9SBarry Smith   c->roworiented = a->roworiented;
1282416022c9SBarry Smith   c->nonew       = a->nonew;
128317ab2063SBarry Smith 
1284416022c9SBarry Smith   if (a->diag) {
12850452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1286416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
128717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1288416022c9SBarry Smith       c->diag[i] = a->diag[i];
128917ab2063SBarry Smith     }
129017ab2063SBarry Smith   }
1291416022c9SBarry Smith   else c->diag        = 0;
1292*754ec7b1SSatish Balay   if( a->inode.size){
1293*754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1294*754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1295*754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1296*754ec7b1SSatish Balay   } else {
1297*754ec7b1SSatish Balay     c->inode.size       = 0;
1298*754ec7b1SSatish Balay     c->inode.node_count = 0;
1299*754ec7b1SSatish Balay   }
1300416022c9SBarry Smith   c->assembled        = 1;
1301416022c9SBarry Smith   c->nz               = a->nz;
1302416022c9SBarry Smith   c->maxnz            = a->maxnz;
1303416022c9SBarry Smith   c->solve_work       = 0;
13044043dd9cSLois Curfman McInnes   c->spptr            = 0;
1305*754ec7b1SSatish Balay   b->inode.node_count = 0;
1306*754ec7b1SSatish Balay   b->inode.size       = 0;
1307*754ec7b1SSatish Balay 
1308416022c9SBarry Smith   *B = C;
130917ab2063SBarry Smith   return 0;
131017ab2063SBarry Smith }
131117ab2063SBarry Smith 
1312416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
131317ab2063SBarry Smith {
1314416022c9SBarry Smith   Mat_SeqAIJ   *a;
1315416022c9SBarry Smith   Mat          B;
131617699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
131717ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
131817ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
131917ab2063SBarry Smith 
132017699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
132117699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
132217ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1323416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
132448d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
132517ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
132617ab2063SBarry Smith 
132717ab2063SBarry Smith   /* read in row lengths */
13280452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1329416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
133017ab2063SBarry Smith 
133117ab2063SBarry Smith   /* create our matrix */
1332416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1333416022c9SBarry Smith   B = *A;
1334416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1335416022c9SBarry Smith   shift = a->indexshift;
133617ab2063SBarry Smith 
133717ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1338416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
133917ab2063SBarry Smith   if (shift) {
134017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1341416022c9SBarry Smith       a->j[i] += 1;
134217ab2063SBarry Smith     }
134317ab2063SBarry Smith   }
134417ab2063SBarry Smith 
134517ab2063SBarry Smith   /* read in nonzero values */
1346416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
134717ab2063SBarry Smith 
134817ab2063SBarry Smith   /* set matrix "i" values */
1349416022c9SBarry Smith   a->i[0] = -shift;
135017ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1351416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1352416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
135317ab2063SBarry Smith   }
13540452661fSBarry Smith   PetscFree(rowlengths);
135517ab2063SBarry Smith 
1356416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1357416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
135817ab2063SBarry Smith   return 0;
135917ab2063SBarry Smith }
136017ab2063SBarry Smith 
136117ab2063SBarry Smith 
136217ab2063SBarry Smith 
1363