xref: /petsc/src/mat/impls/aij/seq/aij.c (revision d5d45c9bad5d66633d6ef030bee5689b038db93c)
117ab2063SBarry Smith #ifndef lint
2*d5d45c9bSBarry Smith static char vcid[] = "$Id: aij.c,v 1.108 1995/11/01 23:17:58 bsmith Exp bsmith $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5*d5d45c9bSBarry Smith /*
6*d5d45c9bSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7*d5d45c9bSBarry Smith   matrix storage format.
8*d5d45c9bSBarry 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;
1817ab2063SBarry Smith   int        ierr, *ia, *ja;
1917ab2063SBarry Smith 
20416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix");
2117ab2063SBarry Smith 
22416022c9SBarry Smith   ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr);
23416022c9SBarry Smith   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
240452661fSBarry Smith   PetscFree(ia); PetscFree(ja);
2517ab2063SBarry Smith   return 0;
2617ab2063SBarry Smith }
2717ab2063SBarry Smith 
28416022c9SBarry Smith #define CHUNKSIZE   10
2917ab2063SBarry Smith 
3017ab2063SBarry Smith /* This version has row oriented v  */
31416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
3217ab2063SBarry Smith {
33416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
34416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
35416022c9SBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen;
36*d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
37416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
3817ab2063SBarry Smith 
3917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
40416022c9SBarry Smith     row  = im[k];
4117ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
42416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
4317ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
4417ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
45416022c9SBarry Smith     low = 0;
4617ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
47416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
48416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
49416022c9SBarry Smith       col = in[l] - shift; value = *v++;
50416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
51416022c9SBarry Smith       while (high-low > 5) {
52416022c9SBarry Smith         t = (low+high)/2;
53416022c9SBarry Smith         if (rp[t] > col) high = t;
54416022c9SBarry Smith         else             low  = t;
5517ab2063SBarry Smith       }
56416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
5717ab2063SBarry Smith         if (rp[i] > col) break;
5817ab2063SBarry Smith         if (rp[i] == col) {
59416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
6017ab2063SBarry Smith           else                  ap[i] = value;
6117ab2063SBarry Smith           goto noinsert;
6217ab2063SBarry Smith         }
6317ab2063SBarry Smith       }
6417ab2063SBarry Smith       if (nonew) goto noinsert;
6517ab2063SBarry Smith       if (nrow >= rmax) {
6617ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
67416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
6817ab2063SBarry Smith         Scalar *new_a;
6917ab2063SBarry Smith 
7017ab2063SBarry Smith         /* malloc new storage space */
71416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
720452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
7317ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
7417ab2063SBarry Smith         new_i   = new_j + new_nz;
7517ab2063SBarry Smith 
7617ab2063SBarry Smith         /* copy over old data into new slots */
7717ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
78416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
79416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
80416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
81416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
8217ab2063SBarry Smith                                                            len*sizeof(int));
83416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
84416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
8517ab2063SBarry Smith                                                            len*sizeof(Scalar));
8617ab2063SBarry Smith         /* free up old matrix storage */
870452661fSBarry Smith         PetscFree(a->a);
880452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
89416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
90416022c9SBarry Smith         a->singlemalloc = 1;
9117ab2063SBarry Smith 
9217ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
93416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
94416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
95416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
9617ab2063SBarry Smith       }
97416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
98416022c9SBarry Smith       /* shift up all the later entries in this row */
99416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
10017ab2063SBarry Smith         rp[ii+1] = rp[ii];
10117ab2063SBarry Smith         ap[ii+1] = ap[ii];
10217ab2063SBarry Smith       }
10317ab2063SBarry Smith       rp[i] = col;
10417ab2063SBarry Smith       ap[i] = value;
10517ab2063SBarry Smith       noinsert:;
106416022c9SBarry Smith       low = i + 1;
10717ab2063SBarry Smith     }
10817ab2063SBarry Smith     ailen[row] = nrow;
10917ab2063SBarry Smith   }
11017ab2063SBarry Smith   return 0;
11117ab2063SBarry Smith }
11217ab2063SBarry Smith 
11317ab2063SBarry Smith #include "draw.h"
11417ab2063SBarry Smith #include "pinclude/pviewer.h"
115416022c9SBarry Smith #include "sysio.h"
11617ab2063SBarry Smith 
117416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
11817ab2063SBarry Smith {
119416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
120416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
12117ab2063SBarry Smith 
122416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1230452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
124416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
125416022c9SBarry Smith   col_lens[1] = a->m;
126416022c9SBarry Smith   col_lens[2] = a->n;
127416022c9SBarry Smith   col_lens[3] = a->nz;
128416022c9SBarry Smith 
129416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
130416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
131416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
13217ab2063SBarry Smith   }
133416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
1340452661fSBarry Smith   PetscFree(col_lens);
135416022c9SBarry Smith 
136416022c9SBarry Smith   /* store column indices (zero start index) */
137416022c9SBarry Smith   if (a->indexshift) {
138416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
13917ab2063SBarry Smith   }
140416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
141416022c9SBarry Smith   if (a->indexshift) {
142416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
14317ab2063SBarry Smith   }
144416022c9SBarry Smith 
145416022c9SBarry Smith   /* store nonzero values */
146416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
14717ab2063SBarry Smith   return 0;
14817ab2063SBarry Smith }
149416022c9SBarry Smith 
150416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
151416022c9SBarry Smith {
152416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
153416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
15417ab2063SBarry Smith   FILE        *fd;
15517ab2063SBarry Smith   char        *outputname;
15617ab2063SBarry Smith 
157416022c9SBarry Smith   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
158416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
159416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
16017ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
16108480c60SBarry Smith     /* no need to print additional information */ ;
16217ab2063SBarry Smith   }
16317ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
16417ab2063SBarry Smith     int nz, nzalloc, mem;
165416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
166416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
16717ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
16817ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
16917ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
17017ab2063SBarry Smith 
17117ab2063SBarry Smith     for (i=0; i<m; i++) {
172416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
17317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
174416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
175416022c9SBarry Smith                    imag(a->a[j]));
17617ab2063SBarry Smith #else
177416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
17817ab2063SBarry Smith #endif
17917ab2063SBarry Smith       }
18017ab2063SBarry Smith     }
18117ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
18217ab2063SBarry Smith   }
18317ab2063SBarry Smith   else {
18417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
18517ab2063SBarry Smith       fprintf(fd,"row %d:",i);
186416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
18717ab2063SBarry Smith #if defined(PETSC_COMPLEX)
188416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
189416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
19017ab2063SBarry Smith         }
19117ab2063SBarry Smith         else {
192416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
19317ab2063SBarry Smith         }
19417ab2063SBarry Smith #else
195416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
19617ab2063SBarry Smith #endif
19717ab2063SBarry Smith       }
19817ab2063SBarry Smith       fprintf(fd,"\n");
19917ab2063SBarry Smith     }
20017ab2063SBarry Smith   }
20117ab2063SBarry Smith   fflush(fd);
202416022c9SBarry Smith   return 0;
203416022c9SBarry Smith }
204416022c9SBarry Smith 
205416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
206416022c9SBarry Smith {
207416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
208cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
209cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
210416022c9SBarry Smith   DrawCtx     draw = (DrawCtx) viewer;
211cddf8d76SBarry Smith   DrawButton  button;
212cddf8d76SBarry Smith 
213416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
214416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
215416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
216416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
217cddf8d76SBarry Smith   color = DRAW_BLUE;
218416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
219cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
220416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
221cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
222cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
223cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
224cddf8d76SBarry Smith #else
225cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
226cddf8d76SBarry Smith #endif
227cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
228cddf8d76SBarry Smith     }
229cddf8d76SBarry Smith   }
230cddf8d76SBarry Smith   color = DRAW_CYAN;
231cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
232cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
233cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
234cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
235cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
236cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
237cddf8d76SBarry Smith     }
238cddf8d76SBarry Smith   }
239cddf8d76SBarry Smith   color = DRAW_RED;
240cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
241cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
242cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
243cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
244cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
245cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
246cddf8d76SBarry Smith #else
247cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
248cddf8d76SBarry Smith #endif
249cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
250416022c9SBarry Smith     }
251416022c9SBarry Smith   }
252416022c9SBarry Smith   DrawFlush(draw);
253cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
254cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
255cddf8d76SBarry Smith 
256cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
257cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
258cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
259cddf8d76SBarry Smith     DrawClear(draw);
260cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
261cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
262cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
263cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
264cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
265cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
266cddf8d76SBarry Smith     w *= scale; h *= scale;
267cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
268cddf8d76SBarry Smith     color = DRAW_BLUE;
269cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
270cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
271cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
272cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
273cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
274cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
275cddf8d76SBarry Smith #else
276cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
277cddf8d76SBarry Smith #endif
278cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
279cddf8d76SBarry Smith       }
280cddf8d76SBarry Smith     }
281cddf8d76SBarry Smith     color = DRAW_CYAN;
282cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
283cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
284cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
285cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
286cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
287cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
288cddf8d76SBarry Smith       }
289cddf8d76SBarry Smith     }
290cddf8d76SBarry Smith     color = DRAW_RED;
291cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
292cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
293cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
294cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
295cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
296cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
297cddf8d76SBarry Smith #else
298cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
299cddf8d76SBarry Smith #endif
300cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
301cddf8d76SBarry Smith       }
302cddf8d76SBarry Smith     }
303cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
304cddf8d76SBarry Smith   }
305416022c9SBarry Smith   return 0;
306416022c9SBarry Smith }
307416022c9SBarry Smith 
308416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
309416022c9SBarry Smith {
310416022c9SBarry Smith   Mat         A = (Mat) obj;
311416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
312416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
313416022c9SBarry Smith 
314416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
315416022c9SBarry Smith   if (!viewer) {
316416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
317416022c9SBarry Smith   }
318416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
319416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
320416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
321416022c9SBarry Smith     }
322416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
323416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
324416022c9SBarry Smith     }
325416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
326416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
327416022c9SBarry Smith     }
328416022c9SBarry Smith   }
329416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
330416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
331416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
33217ab2063SBarry Smith   }
33317ab2063SBarry Smith   return 0;
33417ab2063SBarry Smith }
33517ab2063SBarry Smith 
336416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
33717ab2063SBarry Smith {
338416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
339416022c9SBarry Smith   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
340416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
341416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
34217ab2063SBarry Smith 
34317ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
34417ab2063SBarry Smith 
34517ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
346416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
34717ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
34817ab2063SBarry Smith     if (fshift) {
349416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
35017ab2063SBarry Smith       N = ailen[i];
35117ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
35217ab2063SBarry Smith         ip[j-fshift] = ip[j];
35317ab2063SBarry Smith         ap[j-fshift] = ap[j];
35417ab2063SBarry Smith       }
35517ab2063SBarry Smith     }
35617ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
35717ab2063SBarry Smith   }
35817ab2063SBarry Smith   if (m) {
35917ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
36017ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
36117ab2063SBarry Smith   }
36217ab2063SBarry Smith   /* reset ilen and imax for each row */
36317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
36417ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
36517ab2063SBarry Smith   }
366416022c9SBarry Smith   a->nz = ai[m] + shift;
36717ab2063SBarry Smith 
36817ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
369416022c9SBarry Smith   if (fshift && a->diag) {
3700452661fSBarry Smith     PetscFree(a->diag);
371416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
372416022c9SBarry Smith     a->diag = 0;
37317ab2063SBarry Smith   }
374416022c9SBarry Smith   a->assembled = 1;
37517ab2063SBarry Smith   return 0;
37617ab2063SBarry Smith }
37717ab2063SBarry Smith 
378416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
37917ab2063SBarry Smith {
380416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
381cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
38217ab2063SBarry Smith   return 0;
38317ab2063SBarry Smith }
384416022c9SBarry Smith 
38517ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
38617ab2063SBarry Smith {
387416022c9SBarry Smith   Mat        A  = (Mat) obj;
388416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
389*d5d45c9bSBarry Smith 
39017ab2063SBarry Smith #if defined(PETSC_LOG)
391416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
39217ab2063SBarry Smith #endif
3930452661fSBarry Smith   PetscFree(a->a);
3940452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
3950452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
3960452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
3970452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
3980452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
3990452661fSBarry Smith   PetscFree(a);
400416022c9SBarry Smith   PLogObjectDestroy(A);
4010452661fSBarry Smith   PetscHeaderDestroy(A);
40217ab2063SBarry Smith   return 0;
40317ab2063SBarry Smith }
40417ab2063SBarry Smith 
405416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
40617ab2063SBarry Smith {
40717ab2063SBarry Smith   return 0;
40817ab2063SBarry Smith }
40917ab2063SBarry Smith 
410416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
41117ab2063SBarry Smith {
412416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
413416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
414416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
415416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
416416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
417e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
418e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
419e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
420e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
421df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
422df719601SLois Curfman McInnes   else if (op == COLUMN_ORIENTED)
4234d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");}
424df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
4254d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
426e2f28af5SBarry Smith   else
4274d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
42817ab2063SBarry Smith   return 0;
42917ab2063SBarry Smith }
43017ab2063SBarry Smith 
431416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
43217ab2063SBarry Smith {
433416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
434416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
43517ab2063SBarry Smith   Scalar     *x, zero = 0.0;
43617ab2063SBarry Smith 
437416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
43817ab2063SBarry Smith   VecSet(&zero,v);
43917ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
440416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
441416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
442416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
443416022c9SBarry Smith       if (a->j[j]+shift == i) {
444416022c9SBarry Smith         x[i] = a->a[j];
44517ab2063SBarry Smith         break;
44617ab2063SBarry Smith       }
44717ab2063SBarry Smith     }
44817ab2063SBarry Smith   }
44917ab2063SBarry Smith   return 0;
45017ab2063SBarry Smith }
45117ab2063SBarry Smith 
45217ab2063SBarry Smith /* -------------------------------------------------------*/
45317ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
45417ab2063SBarry Smith /* -------------------------------------------------------*/
455416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
45617ab2063SBarry Smith {
457416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
45817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
459416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
46017ab2063SBarry Smith 
461416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
46217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
463cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
46417ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
46517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
466416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
467416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
468416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
46917ab2063SBarry Smith     alpha = x[i];
47017ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
47117ab2063SBarry Smith   }
472416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
47317ab2063SBarry Smith   return 0;
47417ab2063SBarry Smith }
475*d5d45c9bSBarry Smith 
476416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
47717ab2063SBarry Smith {
478416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
47917ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
480416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
48117ab2063SBarry Smith 
482416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
48317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
48417ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
48517ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
48617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
487416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
488416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
489416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
49017ab2063SBarry Smith     alpha = x[i];
49117ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
49217ab2063SBarry Smith   }
49317ab2063SBarry Smith   return 0;
49417ab2063SBarry Smith }
49517ab2063SBarry Smith 
496416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
49717ab2063SBarry Smith {
498416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
49917ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
500416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
50117ab2063SBarry Smith 
502416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
50317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
50417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
505416022c9SBarry Smith   idx  = a->j;
506416022c9SBarry Smith   v    = a->a;
507416022c9SBarry Smith   ii   = a->i;
50817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
509416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
51017ab2063SBarry Smith     sum  = 0.0;
51117ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
51217ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
513416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
51417ab2063SBarry Smith     y[i] = sum;
51517ab2063SBarry Smith   }
516416022c9SBarry Smith   PLogFlops(2*a->nz - m);
51717ab2063SBarry Smith   return 0;
51817ab2063SBarry Smith }
51917ab2063SBarry Smith 
520416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
52117ab2063SBarry Smith {
522416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
52317ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
524cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
52517ab2063SBarry Smith 
52648d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
52717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
52817ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
529cddf8d76SBarry Smith   idx  = a->j;
530cddf8d76SBarry Smith   v    = a->a;
531cddf8d76SBarry Smith   ii   = a->i;
53217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
533cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
53417ab2063SBarry Smith     sum  = y[i];
535cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
536cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
53717ab2063SBarry Smith     z[i] = sum;
53817ab2063SBarry Smith   }
539416022c9SBarry Smith   PLogFlops(2*a->nz);
54017ab2063SBarry Smith   return 0;
54117ab2063SBarry Smith }
54217ab2063SBarry Smith 
54317ab2063SBarry Smith /*
54417ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
54517ab2063SBarry Smith */
54617ab2063SBarry Smith 
547416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
54817ab2063SBarry Smith {
549416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
550416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
55117ab2063SBarry Smith 
552416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
5530452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
554416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
555416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
556416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
557416022c9SBarry Smith       if (a->j[j]+shift == i) {
55817ab2063SBarry Smith         diag[i] = j - shift;
55917ab2063SBarry Smith         break;
56017ab2063SBarry Smith       }
56117ab2063SBarry Smith     }
56217ab2063SBarry Smith   }
563416022c9SBarry Smith   a->diag = diag;
56417ab2063SBarry Smith   return 0;
56517ab2063SBarry Smith }
56617ab2063SBarry Smith 
567416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
56817ab2063SBarry Smith                            double fshift,int its,Vec xx)
56917ab2063SBarry Smith {
570416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
571416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
572*d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
57317ab2063SBarry Smith 
57417ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
575416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
576416022c9SBarry Smith   diag = a->diag;
577416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
57817ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
57917ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
58017ab2063SBarry Smith     bs = b + shift;
58117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
582416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
583416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
584416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
585416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
58617ab2063SBarry Smith         sum  = b[i]*d/omega;
58717ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
58817ab2063SBarry Smith         x[i] = sum;
58917ab2063SBarry Smith     }
59017ab2063SBarry Smith     return 0;
59117ab2063SBarry Smith   }
59217ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
59317ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
59417ab2063SBarry Smith   }
595416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
59617ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
59717ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
59817ab2063SBarry Smith 
59917ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
60017ab2063SBarry Smith 
60117ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
60217ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
60317ab2063SBarry Smith     is the relaxation factor.
60417ab2063SBarry Smith     */
6050452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
60617ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
60717ab2063SBarry Smith 
60817ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
60917ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
610416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
611416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
612416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
613416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
61417ab2063SBarry Smith       sum  = b[i];
61517ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
61617ab2063SBarry Smith       x[i] = omega*(sum/d);
61717ab2063SBarry Smith     }
61817ab2063SBarry Smith 
61917ab2063SBarry Smith     /*  t = b - (2*E - D)x */
620416022c9SBarry Smith     v = a->a;
62117ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
62217ab2063SBarry Smith 
62317ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
624416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
625416022c9SBarry Smith     diag = a->diag;
62617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
627416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
628416022c9SBarry Smith       n    = diag[i] - a->i[i];
629416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
630416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
63117ab2063SBarry Smith       sum  = t[i];
63217ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
63317ab2063SBarry Smith       t[i] = omega*(sum/d);
63417ab2063SBarry Smith     }
63517ab2063SBarry Smith 
63617ab2063SBarry Smith     /*  x = x + t */
63717ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
6380452661fSBarry Smith     PetscFree(t);
63917ab2063SBarry Smith     return 0;
64017ab2063SBarry Smith   }
64117ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
64217ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
64317ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
644416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
645416022c9SBarry Smith         n    = diag[i] - a->i[i];
646416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
647416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
64817ab2063SBarry Smith         sum  = b[i];
64917ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
65017ab2063SBarry Smith         x[i] = omega*(sum/d);
65117ab2063SBarry Smith       }
65217ab2063SBarry Smith       xb = x;
65317ab2063SBarry Smith     }
65417ab2063SBarry Smith     else xb = b;
65517ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
65617ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
65717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
658416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
65917ab2063SBarry Smith       }
66017ab2063SBarry Smith     }
66117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
66217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
663416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
664416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
665416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
666416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
66717ab2063SBarry Smith         sum  = xb[i];
66817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
66917ab2063SBarry Smith         x[i] = omega*(sum/d);
67017ab2063SBarry Smith       }
67117ab2063SBarry Smith     }
67217ab2063SBarry Smith     its--;
67317ab2063SBarry Smith   }
67417ab2063SBarry Smith   while (its--) {
67517ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
67617ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
677416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
678416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
679416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
680416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
68117ab2063SBarry Smith         sum  = b[i];
68217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
68317ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
68417ab2063SBarry Smith       }
68517ab2063SBarry Smith     }
68617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
68717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
688416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
689416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
690416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
691416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
69217ab2063SBarry Smith         sum  = b[i];
69317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
69417ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
69517ab2063SBarry Smith       }
69617ab2063SBarry Smith     }
69717ab2063SBarry Smith   }
69817ab2063SBarry Smith   return 0;
69917ab2063SBarry Smith }
70017ab2063SBarry Smith 
701*d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
70217ab2063SBarry Smith {
703416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
704416022c9SBarry Smith   *nz      = a->nz;
705416022c9SBarry Smith   *nzalloc = a->maxnz;
706416022c9SBarry Smith   *mem     = (int)A->mem;
70717ab2063SBarry Smith   return 0;
70817ab2063SBarry Smith }
70917ab2063SBarry Smith 
71017ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
71117ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
71217ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
71317ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
71417ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
71517ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
71617ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
71717ab2063SBarry Smith 
71817ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
71917ab2063SBarry Smith {
720416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
721416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
72217ab2063SBarry Smith 
72317ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
72417ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
72517ab2063SBarry Smith   if (diag) {
72617ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
727416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
728416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
729416022c9SBarry Smith         a->ilen[rows[i]] = 1;
730416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
731416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
73217ab2063SBarry Smith       }
73317ab2063SBarry Smith       else {
73417ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
73517ab2063SBarry Smith         CHKERRQ(ierr);
73617ab2063SBarry Smith       }
73717ab2063SBarry Smith     }
73817ab2063SBarry Smith   }
73917ab2063SBarry Smith   else {
74017ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
741416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
742416022c9SBarry Smith       a->ilen[rows[i]] = 0;
74317ab2063SBarry Smith     }
74417ab2063SBarry Smith   }
74517ab2063SBarry Smith   ISRestoreIndices(is,&rows);
74617ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
74717ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
74817ab2063SBarry Smith   return 0;
74917ab2063SBarry Smith }
75017ab2063SBarry Smith 
751416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
75217ab2063SBarry Smith {
753416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
754416022c9SBarry Smith   *m = a->m; *n = a->n;
75517ab2063SBarry Smith   return 0;
75617ab2063SBarry Smith }
75717ab2063SBarry Smith 
758416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
75917ab2063SBarry Smith {
760416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
761416022c9SBarry Smith   *m = 0; *n = a->m;
76217ab2063SBarry Smith   return 0;
76317ab2063SBarry Smith }
764416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
76517ab2063SBarry Smith {
766416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
767416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
76817ab2063SBarry Smith 
769416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
77017ab2063SBarry Smith 
771416022c9SBarry Smith   if (!a->assembled) {
772416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
773416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
77417ab2063SBarry Smith   }
775416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
776416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
77717ab2063SBarry Smith   if (idx) {
77817ab2063SBarry Smith     if (*nz) {
779416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
7800452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
78117ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
78217ab2063SBarry Smith     }
78317ab2063SBarry Smith     else *idx = 0;
78417ab2063SBarry Smith   }
78517ab2063SBarry Smith   return 0;
78617ab2063SBarry Smith }
78717ab2063SBarry Smith 
788416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
78917ab2063SBarry Smith {
7900452661fSBarry Smith   if (idx) {if (*idx) PetscFree(*idx);}
79117ab2063SBarry Smith   return 0;
79217ab2063SBarry Smith }
79317ab2063SBarry Smith 
794cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
79517ab2063SBarry Smith {
796416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
797416022c9SBarry Smith   Scalar     *v = a->a;
79817ab2063SBarry Smith   double     sum = 0.0;
799416022c9SBarry Smith   int        i, j,shift = a->indexshift;
80017ab2063SBarry Smith 
801416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
80217ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
803416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
80417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
80517ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
80617ab2063SBarry Smith #else
80717ab2063SBarry Smith       sum += (*v)*(*v); v++;
80817ab2063SBarry Smith #endif
80917ab2063SBarry Smith     }
81017ab2063SBarry Smith     *norm = sqrt(sum);
81117ab2063SBarry Smith   }
81217ab2063SBarry Smith   else if (type == NORM_1) {
81317ab2063SBarry Smith     double *tmp;
814416022c9SBarry Smith     int    *jj = a->j;
8150452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
816cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
81717ab2063SBarry Smith     *norm = 0.0;
818416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
819cddf8d76SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v++);
82017ab2063SBarry Smith     }
821416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
82217ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
82317ab2063SBarry Smith     }
8240452661fSBarry Smith     PetscFree(tmp);
82517ab2063SBarry Smith   }
82617ab2063SBarry Smith   else if (type == NORM_INFINITY) {
82717ab2063SBarry Smith     *norm = 0.0;
828416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
829416022c9SBarry Smith       v = a->a + a->i[j] + shift;
83017ab2063SBarry Smith       sum = 0.0;
831416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
832cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
83317ab2063SBarry Smith       }
83417ab2063SBarry Smith       if (sum > *norm) *norm = sum;
83517ab2063SBarry Smith     }
83617ab2063SBarry Smith   }
83717ab2063SBarry Smith   else {
83848d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
83917ab2063SBarry Smith   }
84017ab2063SBarry Smith   return 0;
84117ab2063SBarry Smith }
84217ab2063SBarry Smith 
843416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
84417ab2063SBarry Smith {
845416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
846416022c9SBarry Smith   Mat        C;
847416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
848416022c9SBarry Smith   int        shift = a->indexshift;
849*d5d45c9bSBarry Smith   Scalar     *array = a->a;
85017ab2063SBarry Smith 
851416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
8520452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
853cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
85417ab2063SBarry Smith   if (shift) {
85517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
85617ab2063SBarry Smith   }
85717ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
858416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
8590452661fSBarry Smith   PetscFree(col);
86017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
86117ab2063SBarry Smith     len = ai[i+1]-ai[i];
862416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
86317ab2063SBarry Smith     array += len; aj += len;
86417ab2063SBarry Smith   }
86517ab2063SBarry Smith   if (shift) {
86617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
86717ab2063SBarry Smith   }
86817ab2063SBarry Smith 
869416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
870416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
87117ab2063SBarry Smith 
872416022c9SBarry Smith   if (B) {
873416022c9SBarry Smith     *B = C;
87417ab2063SBarry Smith   } else {
875416022c9SBarry Smith     /* This isn't really an in-place transpose */
8760452661fSBarry Smith     PetscFree(a->a);
8770452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
8780452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
8790452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
8800452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
8810452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
8820452661fSBarry Smith     PetscFree(a);
883416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
8840452661fSBarry Smith     PetscHeaderDestroy(C);
88517ab2063SBarry Smith   }
88617ab2063SBarry Smith   return 0;
88717ab2063SBarry Smith }
88817ab2063SBarry Smith 
889416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
89017ab2063SBarry Smith {
891416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
89217ab2063SBarry Smith   Scalar     *l,*r,x,*v;
893*d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
89417ab2063SBarry Smith 
89548d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
89617ab2063SBarry Smith   if (ll) {
89717ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
898416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
899416022c9SBarry Smith     v = a->a;
90017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
90117ab2063SBarry Smith       x = l[i];
902416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
90317ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
90417ab2063SBarry Smith     }
90517ab2063SBarry Smith   }
90617ab2063SBarry Smith   if (rr) {
90717ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
908416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
909416022c9SBarry Smith     v = a->a; jj = a->j;
91017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
91117ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
91217ab2063SBarry Smith     }
91317ab2063SBarry Smith   }
91417ab2063SBarry Smith   return 0;
91517ab2063SBarry Smith }
91617ab2063SBarry Smith 
917cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
91817ab2063SBarry Smith {
919db02288aSLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c;
92002834360SBarry Smith   int        nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
92102834360SBarry Smith   int        *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
922db02288aSLois Curfman McInnes   int        first,step,*starts,*j_new,*i_new;
923db02288aSLois Curfman McInnes   Scalar     *vwork,*a_new;
924416022c9SBarry Smith   Mat        C;
92517ab2063SBarry Smith 
926416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
92717ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
92817ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
92917ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
93017ab2063SBarry Smith 
93102834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
93202834360SBarry Smith     /* special case of contiguous rows */
9330452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
93402834360SBarry Smith     starts = lens + ncols;
93502834360SBarry Smith     /* loop over new rows determining lens and starting points */
93602834360SBarry Smith     for (i=0; i<nrows; i++) {
93702834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
93802834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
93902834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
94002834360SBarry Smith         if (a->j[k] >= first) {
94102834360SBarry Smith           starts[i] = k;
94202834360SBarry Smith           break;
94302834360SBarry Smith 	}
94402834360SBarry Smith       }
94502834360SBarry Smith       lens[i] = 0;
94602834360SBarry Smith       while (k < kend) {
94702834360SBarry Smith         if (a->j[k++] >= first+ncols) break;
94802834360SBarry Smith         lens[i]++;
94902834360SBarry Smith       }
95002834360SBarry Smith     }
95102834360SBarry Smith     /* create submatrix */
952cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
95308480c60SBarry Smith       int n_cols,n_rows;
95408480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
95508480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
95608480c60SBarry Smith       MatZeroEntries(*B);
95708480c60SBarry Smith       C = *B;
95808480c60SBarry Smith     }
95908480c60SBarry Smith     else {
96002834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
96108480c60SBarry Smith     }
962db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
963db02288aSLois Curfman McInnes 
96402834360SBarry Smith     /* loop over rows inserting into submatrix */
965db02288aSLois Curfman McInnes     a_new    = c->a;
966db02288aSLois Curfman McInnes     j_new    = c->j;
967db02288aSLois Curfman McInnes     i_new    = c->i;
968db02288aSLois Curfman McInnes     i_new[0] = -shift;
96902834360SBarry Smith     for (i=0; i<nrows; i++) {
97002834360SBarry Smith       for ( k=0; k<lens[i]; k++ ) {
971db02288aSLois Curfman McInnes         *j_new++ = a->j[k+starts[i]] - first;
97202834360SBarry Smith       }
973db02288aSLois Curfman McInnes       PetscMemcpy(a_new,a->a + starts[i],lens[i]*sizeof(Scalar));
974db02288aSLois Curfman McInnes       a_new      += lens[i];
975db02288aSLois Curfman McInnes       i_new[i+1]  = i_new[i] + lens[i];
9761987afe7SBarry Smith       c->ilen[i]  = lens[i];
97702834360SBarry Smith     }
9780452661fSBarry Smith     PetscFree(lens);
97902834360SBarry Smith   }
98002834360SBarry Smith   else {
98102834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
9820452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
98302834360SBarry Smith     ssmap = smap + shift;
9840452661fSBarry Smith     cwork = (int *) PetscMalloc((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork);
98502834360SBarry Smith     lens  = cwork + ncols;
9860452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
987cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
98817ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
98902834360SBarry Smith     /* determine lens of each row */
99002834360SBarry Smith     for (i=0; i<nrows; i++) {
99102834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
99202834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
99302834360SBarry Smith       lens[i] = 0;
99402834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
99502834360SBarry Smith         if (ssmap[a->j[k]]) {
99602834360SBarry Smith           lens[i]++;
99702834360SBarry Smith         }
99802834360SBarry Smith       }
99902834360SBarry Smith     }
100017ab2063SBarry Smith     /* Create and fill new matrix */
100108480c60SBarry Smith     if (MatValidMatrix(*B)) {
100208480c60SBarry Smith       int n_cols,n_rows;
100308480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
100408480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
100508480c60SBarry Smith       MatZeroEntries(*B);
100608480c60SBarry Smith       C = *B;
100708480c60SBarry Smith     }
100808480c60SBarry Smith     else {
100902834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
101008480c60SBarry Smith     }
101117ab2063SBarry Smith     for (i=0; i<nrows; i++) {
101217ab2063SBarry Smith       nznew  = 0;
1013416022c9SBarry Smith       kstart = a->i[irow[i]]+shift;
1014416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
101517ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
101602834360SBarry Smith         if (ssmap[a->j[k]]) {
101702834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1018416022c9SBarry Smith           vwork[nznew++] = a->a[k];
101917ab2063SBarry Smith         }
102017ab2063SBarry Smith       }
1021416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
102217ab2063SBarry Smith     }
102302834360SBarry Smith     /* Free work space */
102402834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10250452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
102602834360SBarry Smith   }
1027416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1028416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
102917ab2063SBarry Smith 
103017ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1031416022c9SBarry Smith   *B = C;
103217ab2063SBarry Smith   return 0;
103317ab2063SBarry Smith }
103417ab2063SBarry Smith 
1035a871dcd8SBarry Smith /*
103663b91edcSBarry Smith      note: This can only work for identity for row and col. It would
103763b91edcSBarry Smith    be good to check this and otherwise generate an error.
1038a871dcd8SBarry Smith */
103963b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1040a871dcd8SBarry Smith {
104163b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
104208480c60SBarry Smith   int        ierr;
104363b91edcSBarry Smith   Mat        outA;
104463b91edcSBarry Smith 
1045a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1046a871dcd8SBarry Smith 
104763b91edcSBarry Smith   outA          = inA;
104863b91edcSBarry Smith   inA->factor   = FACTOR_LU;
104963b91edcSBarry Smith   a->row        = row;
105063b91edcSBarry Smith   a->col        = col;
105163b91edcSBarry Smith 
10520452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
105363b91edcSBarry Smith 
105408480c60SBarry Smith   if (!a->diag) {
105508480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
105663b91edcSBarry Smith   }
105763b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1058a871dcd8SBarry Smith   return 0;
1059a871dcd8SBarry Smith }
1060a871dcd8SBarry Smith 
1061cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1062cddf8d76SBarry Smith                                     Mat **B)
1063cddf8d76SBarry Smith {
1064cddf8d76SBarry Smith   int ierr,i;
1065cddf8d76SBarry Smith 
1066cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
10670452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1068cddf8d76SBarry Smith   }
1069cddf8d76SBarry Smith 
1070cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1071cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1072cddf8d76SBarry Smith   }
1073cddf8d76SBarry Smith   return 0;
1074cddf8d76SBarry Smith }
1075cddf8d76SBarry Smith 
107617ab2063SBarry Smith /* -------------------------------------------------------------------*/
107717ab2063SBarry Smith 
107817ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
107917ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1080416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1081416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
108217ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
108317ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
108417ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
108517ab2063SBarry Smith        MatRelax_SeqAIJ,
108617ab2063SBarry Smith        MatTranspose_SeqAIJ,
108717ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
108817ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
108917ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
109017ab2063SBarry Smith        MatCompress_SeqAIJ,
109117ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
109217ab2063SBarry Smith        MatGetReordering_SeqAIJ,
109317ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
109417ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
109517ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
109617ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1097416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
1098a871dcd8SBarry Smith        MatCopyPrivate_SeqAIJ,0,0,
1099cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
1100cddf8d76SBarry Smith        MatGetSubMatrices_SeqAIJ};
110117ab2063SBarry Smith 
110217ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
110317ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
110417ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
110517ab2063SBarry Smith 
110617ab2063SBarry Smith /*@C
110717ab2063SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
110817ab2063SBarry Smith    (the default uniprocessor PETSc format).
110917ab2063SBarry Smith 
111017ab2063SBarry Smith    Input Parameters:
111117ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
111217ab2063SBarry Smith .  m - number of rows
111317ab2063SBarry Smith .  n - number of columns
111417ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
111517ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
111617ab2063SBarry Smith 
111717ab2063SBarry Smith    Output Parameter:
1118416022c9SBarry Smith .  A - the matrix
111917ab2063SBarry Smith 
112017ab2063SBarry Smith    Notes:
112117ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
112217ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
11230002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
11240002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
112517ab2063SBarry Smith 
112617ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
112717ab2063SBarry Smith    Set both nz and nnz to zero for PETSc to control dynamic memory
112817ab2063SBarry Smith    allocation.
112917ab2063SBarry Smith 
113017ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse
113117ab2063SBarry Smith 
113217ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
113317ab2063SBarry Smith @*/
1134416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
113517ab2063SBarry Smith {
1136416022c9SBarry Smith   Mat        B;
1137416022c9SBarry Smith   Mat_SeqAIJ *b;
113817ab2063SBarry Smith   int        i,len,ierr;
1139*d5d45c9bSBarry Smith 
1140416022c9SBarry Smith   *A      = 0;
11410452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1142416022c9SBarry Smith   PLogObjectCreate(B);
11430452661fSBarry Smith   B->data               = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1144416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1145416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1146416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1147416022c9SBarry Smith   B->factor           = 0;
1148416022c9SBarry Smith   B->lupivotthreshold = 1.0;
1149416022c9SBarry Smith   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1150416022c9SBarry Smith   b->row              = 0;
1151416022c9SBarry Smith   b->col              = 0;
1152416022c9SBarry Smith   b->indexshift       = 0;
1153416022c9SBarry Smith   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
115417ab2063SBarry Smith 
1155416022c9SBarry Smith   b->m       = m;
1156416022c9SBarry Smith   b->n       = n;
11570452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
115817ab2063SBarry Smith   if (!nnz) {
115917ab2063SBarry Smith     if (nz <= 0) nz = 1;
1160416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
116117ab2063SBarry Smith     nz = nz*m;
116217ab2063SBarry Smith   }
116317ab2063SBarry Smith   else {
116417ab2063SBarry Smith     nz = 0;
1165416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
116617ab2063SBarry Smith   }
116717ab2063SBarry Smith 
116817ab2063SBarry Smith   /* allocate the matrix space */
1169416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
11700452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1171416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1172cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1173416022c9SBarry Smith   b->i  = b->j + nz;
1174416022c9SBarry Smith   b->singlemalloc = 1;
117517ab2063SBarry Smith 
1176416022c9SBarry Smith   b->i[0] = -b->indexshift;
117717ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1178416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
117917ab2063SBarry Smith   }
118017ab2063SBarry Smith 
1181416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
11820452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1183416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1184416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
118517ab2063SBarry Smith 
1186416022c9SBarry Smith   b->nz          = 0;
1187416022c9SBarry Smith   b->maxnz       = nz;
1188416022c9SBarry Smith   b->sorted      = 0;
1189416022c9SBarry Smith   b->roworiented = 1;
1190416022c9SBarry Smith   b->nonew       = 0;
1191416022c9SBarry Smith   b->diag        = 0;
1192416022c9SBarry Smith   b->assembled   = 0;
1193416022c9SBarry Smith   b->solve_work  = 0;
119471bd300dSLois Curfman McInnes   b->spptr       = 0;
119517ab2063SBarry Smith 
1196416022c9SBarry Smith   *A = B;
119717ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_superlu")) {
1198416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
119917ab2063SBarry Smith   }
120017ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_essl")) {
1201416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
120217ab2063SBarry Smith   }
120317ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_dxml")) {
1204416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1205416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
120617ab2063SBarry Smith   }
120717ab2063SBarry Smith 
120817ab2063SBarry Smith   return 0;
120917ab2063SBarry Smith }
121017ab2063SBarry Smith 
121108480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues)
121217ab2063SBarry Smith {
1213416022c9SBarry Smith   Mat        C;
1214416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
121508480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
121617ab2063SBarry Smith 
12174043dd9cSLois Curfman McInnes   *B = 0;
1218416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
12190452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1220416022c9SBarry Smith   PLogObjectCreate(C);
12210452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1222416022c9SBarry Smith   PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps));
1223416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1224416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1225416022c9SBarry Smith   C->factor     = A->factor;
1226416022c9SBarry Smith   c->row        = 0;
1227416022c9SBarry Smith   c->col        = 0;
1228416022c9SBarry Smith   c->indexshift = shift;
122917ab2063SBarry Smith 
1230416022c9SBarry Smith   c->m          = a->m;
1231416022c9SBarry Smith   c->n          = a->n;
123217ab2063SBarry Smith 
12330452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
12340452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
123517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1236416022c9SBarry Smith     c->imax[i] = a->imax[i];
1237416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
123817ab2063SBarry Smith   }
123917ab2063SBarry Smith 
124017ab2063SBarry Smith   /* allocate the matrix space */
1241416022c9SBarry Smith   c->singlemalloc = 1;
1242416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
12430452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1244416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1245416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1246416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
124717ab2063SBarry Smith   if (m > 0) {
1248416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
124908480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1250416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
125117ab2063SBarry Smith     }
125208480c60SBarry Smith   }
125317ab2063SBarry Smith 
1254416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1255416022c9SBarry Smith   c->sorted      = a->sorted;
1256416022c9SBarry Smith   c->roworiented = a->roworiented;
1257416022c9SBarry Smith   c->nonew       = a->nonew;
125817ab2063SBarry Smith 
1259416022c9SBarry Smith   if (a->diag) {
12600452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1261416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
126217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1263416022c9SBarry Smith       c->diag[i] = a->diag[i];
126417ab2063SBarry Smith     }
126517ab2063SBarry Smith   }
1266416022c9SBarry Smith   else c->diag        = 0;
1267416022c9SBarry Smith   c->assembled        = 1;
1268416022c9SBarry Smith   c->nz               = a->nz;
1269416022c9SBarry Smith   c->maxnz            = a->maxnz;
1270416022c9SBarry Smith   c->solve_work       = 0;
12714043dd9cSLois Curfman McInnes   c->spptr            = 0;
1272416022c9SBarry Smith   *B = C;
127317ab2063SBarry Smith   return 0;
127417ab2063SBarry Smith }
127517ab2063SBarry Smith 
1276416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
127717ab2063SBarry Smith {
1278416022c9SBarry Smith   Mat_SeqAIJ   *a;
1279416022c9SBarry Smith   Mat          B;
128017699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
128117ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
128217ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
128317ab2063SBarry Smith 
128417699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
128517699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
128617ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1287416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
128848d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
128917ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
129017ab2063SBarry Smith 
129117ab2063SBarry Smith   /* read in row lengths */
12920452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1293416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
129417ab2063SBarry Smith 
129517ab2063SBarry Smith   /* create our matrix */
1296416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1297416022c9SBarry Smith   B = *A;
1298416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1299416022c9SBarry Smith   shift = a->indexshift;
130017ab2063SBarry Smith 
130117ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1302416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
130317ab2063SBarry Smith   if (shift) {
130417ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1305416022c9SBarry Smith       a->j[i] += 1;
130617ab2063SBarry Smith     }
130717ab2063SBarry Smith   }
130817ab2063SBarry Smith 
130917ab2063SBarry Smith   /* read in nonzero values */
1310416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
131117ab2063SBarry Smith 
131217ab2063SBarry Smith   /* set matrix "i" values */
1313416022c9SBarry Smith   a->i[0] = -shift;
131417ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1315416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1316416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
131717ab2063SBarry Smith   }
13180452661fSBarry Smith   PetscFree(rowlengths);
131917ab2063SBarry Smith 
1320416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1321416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
132217ab2063SBarry Smith   return 0;
132317ab2063SBarry Smith }
132417ab2063SBarry Smith 
132517ab2063SBarry Smith 
132617ab2063SBarry Smith 
1327