xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 13a67aface5a148d688de3dc83f9ec17fbd665a2)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: aij.c,v 1.167 1996/04/05 16:56:58 balay Exp balay $";
4 #endif
5 
6 /*
7     Defines the basic matrix operations for the AIJ (compressed row)
8   matrix storage format.
9 */
10 #include "aij.h"
11 #include "vec/vecimpl.h"
12 #include "inline/spops.h"
13 #include "petsc.h"
14 #include "inline/bitarray.h"
15 
16 extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
17 
18 static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
19 {
20   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
21   int        ierr, *ia, *ja,n,*idx,i,oshift,ishift;
22 
23   /*
24      this is tacky: In the future when we have written special factorization
25      and solve routines for the identity permutation we should use a
26      stride index set instead of the general one.
27   */
28   if (type  == ORDER_NATURAL) {
29     n = a->n;
30     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
31     for ( i=0; i<n; i++ ) idx[i] = i;
32     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
33     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
34     PetscFree(idx);
35     ISSetPermutation(*rperm);
36     ISSetPermutation(*cperm);
37     ISSetIdentity(*rperm);
38     ISSetIdentity(*cperm);
39     return 0;
40   }
41 
42   MatReorderingRegisterAll();
43   ishift = a->indexshift;
44   oshift = -MatReorderingIndexShift[(int)type];
45   if (MatReorderingRequiresSymmetric[(int)type]) {
46     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja);
47     CHKERRQ(ierr);
48     ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
49     PetscFree(ia); PetscFree(ja);
50   } else {
51     if (ishift == oshift) {
52       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
53     }
54     else if (ishift == -1) {
55       /* temporarily subtract 1 from i and j indices */
56       int nz = a->i[a->n] - 1;
57       for ( i=0; i<nz; i++ ) a->j[i]--;
58       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
59       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
60       for ( i=0; i<nz; i++ ) a->j[i]++;
61       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
62     } else {
63       /* temporarily add 1 to i and j indices */
64       int nz = a->i[a->n] - 1;
65       for ( i=0; i<nz; i++ ) a->j[i]++;
66       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
67       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
68       for ( i=0; i<nz; i++ ) a->j[i]--;
69       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
70     }
71   }
72   return 0;
73 }
74 
75 #define CHUNKSIZE   15
76 
77 /* This version has row oriented v  */
78 static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
79 {
80   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
81   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
82   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
83   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
84   Scalar     *ap,value, *aa = a->a;
85 
86   for ( k=0; k<m; k++ ) { /* loop over added rows */
87     row  = im[k];
88     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
89     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
90     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
91     rmax = imax[row]; nrow = ailen[row];
92     low = 0;
93     for ( l=0; l<n; l++ ) { /* loop over added columns */
94       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
95       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
96       col = in[l] - shift;
97       if (roworiented) {
98         value = *v++;
99       }
100       else {
101         value = v[k + l*m];
102       }
103       if (!sorted) low = 0; high = nrow;
104       while (high-low > 5) {
105         t = (low+high)/2;
106         if (rp[t] > col) high = t;
107         else             low  = t;
108       }
109       for ( i=low; i<high; i++ ) {
110         if (rp[i] > col) break;
111         if (rp[i] == col) {
112           if (is == ADD_VALUES) ap[i] += value;
113           else                  ap[i] = value;
114           goto noinsert;
115         }
116       }
117       if (nonew) goto noinsert;
118       if (nrow >= rmax) {
119         /* there is no extra room in row, therefore enlarge */
120         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
121         Scalar *new_a;
122 
123         /* malloc new storage space */
124         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
125         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
126         new_j   = (int *) (new_a + new_nz);
127         new_i   = new_j + new_nz;
128 
129         /* copy over old data into new slots */
130         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
131         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
132         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
133         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
134         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
135                                                            len*sizeof(int));
136         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
137         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
138                                                            len*sizeof(Scalar));
139         /* free up old matrix storage */
140         PetscFree(a->a);
141         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
142         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
143         a->singlemalloc = 1;
144 
145         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
146         rmax = imax[row] = imax[row] + CHUNKSIZE;
147         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
148         a->maxnz += CHUNKSIZE;
149         a->reallocs++;
150       }
151       N = nrow++ - 1; a->nz++;
152       /* shift up all the later entries in this row */
153       for ( ii=N; ii>=i; ii-- ) {
154         rp[ii+1] = rp[ii];
155         ap[ii+1] = ap[ii];
156       }
157       rp[i] = col;
158       ap[i] = value;
159       noinsert:;
160       low = i + 1;
161     }
162     ailen[row] = nrow;
163   }
164   return 0;
165 }
166 
167 static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
168 {
169   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
170   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
171   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
172   Scalar     *ap, *aa = a->a, zero = 0.0;
173 
174   for ( k=0; k<m; k++ ) { /* loop over rows */
175     row  = im[k];
176     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
177     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
178     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
179     nrow = ailen[row];
180     for ( l=0; l<n; l++ ) { /* loop over columns */
181       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
182       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
183       col = in[l] - shift;
184       high = nrow; low = 0; /* assume unsorted */
185       while (high-low > 5) {
186         t = (low+high)/2;
187         if (rp[t] > col) high = t;
188         else             low  = t;
189       }
190       for ( i=low; i<high; i++ ) {
191         if (rp[i] > col) break;
192         if (rp[i] == col) {
193           *v++ = ap[i];
194           goto finished;
195         }
196       }
197       *v++ = zero;
198       finished:;
199     }
200   }
201   return 0;
202 }
203 
204 #include "draw.h"
205 #include "pinclude/pviewer.h"
206 #include "sys.h"
207 
208 static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
209 {
210   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
211   int        i, fd, *col_lens, ierr;
212 
213   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
214   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
215   col_lens[0] = MAT_COOKIE;
216   col_lens[1] = a->m;
217   col_lens[2] = a->n;
218   col_lens[3] = a->nz;
219 
220   /* store lengths of each row and write (including header) to file */
221   for ( i=0; i<a->m; i++ ) {
222     col_lens[4+i] = a->i[i+1] - a->i[i];
223   }
224   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
225   PetscFree(col_lens);
226 
227   /* store column indices (zero start index) */
228   if (a->indexshift) {
229     for ( i=0; i<a->nz; i++ ) a->j[i]--;
230   }
231   ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr);
232   if (a->indexshift) {
233     for ( i=0; i<a->nz; i++ ) a->j[i]++;
234   }
235 
236   /* store nonzero values */
237   ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
238   return 0;
239 }
240 
241 static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
242 {
243   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
244   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg;
245   FILE        *fd;
246   char        *outputname;
247 
248   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
249   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
250   ierr = ViewerGetFormat(viewer,&format);
251   if (format == ASCII_FORMAT_INFO) {
252     return 0;
253   }
254   else if (format == ASCII_FORMAT_INFO_DETAILED) {
255     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr);
256     if (flg) fprintf(fd,"  not using I-node routines\n");
257     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
258         a->inode.node_count,a->inode.limit);
259   }
260   else if (format == ASCII_FORMAT_MATLAB) {
261     int nz, nzalloc, mem;
262     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
263     fprintf(fd,"%% Size = %d %d \n",m,a->n);
264     fprintf(fd,"%% Nonzeros = %d \n",nz);
265     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
266     fprintf(fd,"zzz = [\n");
267 
268     for (i=0; i<m; i++) {
269       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
270 #if defined(PETSC_COMPLEX)
271         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]),
272                    imag(a->a[j]));
273 #else
274         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
275 #endif
276       }
277     }
278     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
279   }
280   else {
281     for ( i=0; i<m; i++ ) {
282       fprintf(fd,"row %d:",i);
283       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
284 #if defined(PETSC_COMPLEX)
285         if (imag(a->a[j]) != 0.0) {
286           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
287         }
288         else {
289           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
290         }
291 #else
292         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
293 #endif
294       }
295       fprintf(fd,"\n");
296     }
297   }
298   fflush(fd);
299   return 0;
300 }
301 
302 static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
303 {
304   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
305   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
306   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
307   Draw        draw;
308   DrawButton  button;
309   PetscTruth  isnull;
310 
311   ViewerDrawGetDraw(viewer,&draw);
312   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
313 
314   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
315   xr += w;    yr += h;  xl = -w;     yl = -h;
316   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
317   /* loop over matrix elements drawing boxes */
318   color = DRAW_BLUE;
319   for ( i=0; i<m; i++ ) {
320     y_l = m - i - 1.0; y_r = y_l + 1.0;
321     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
322       x_l = a->j[j] + shift; x_r = x_l + 1.0;
323 #if defined(PETSC_COMPLEX)
324       if (real(a->a[j]) >=  0.) continue;
325 #else
326       if (a->a[j] >=  0.) continue;
327 #endif
328       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
329     }
330   }
331   color = DRAW_CYAN;
332   for ( i=0; i<m; i++ ) {
333     y_l = m - i - 1.0; y_r = y_l + 1.0;
334     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
335       x_l = a->j[j] + shift; x_r = x_l + 1.0;
336       if (a->a[j] !=  0.) continue;
337       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
338     }
339   }
340   color = DRAW_RED;
341   for ( i=0; i<m; i++ ) {
342     y_l = m - i - 1.0; y_r = y_l + 1.0;
343     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
344       x_l = a->j[j] + shift; x_r = x_l + 1.0;
345 #if defined(PETSC_COMPLEX)
346       if (real(a->a[j]) <=  0.) continue;
347 #else
348       if (a->a[j] <=  0.) continue;
349 #endif
350       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
351     }
352   }
353   DrawFlush(draw);
354   DrawGetPause(draw,&pause);
355   if (pause >= 0) { PetscSleep(pause); return 0;}
356 
357   /* allow the matrix to zoom or shrink */
358   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
359   while (button != BUTTON_RIGHT) {
360     DrawClear(draw);
361     if (button == BUTTON_LEFT) scale = .5;
362     else if (button == BUTTON_CENTER) scale = 2.;
363     xl = scale*(xl + w - xc) + xc - w*scale;
364     xr = scale*(xr - w - xc) + xc + w*scale;
365     yl = scale*(yl + h - yc) + yc - h*scale;
366     yr = scale*(yr - h - yc) + yc + h*scale;
367     w *= scale; h *= scale;
368     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
369     color = DRAW_BLUE;
370     for ( i=0; i<m; i++ ) {
371       y_l = m - i - 1.0; y_r = y_l + 1.0;
372       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
373         x_l = a->j[j] + shift; x_r = x_l + 1.0;
374 #if defined(PETSC_COMPLEX)
375         if (real(a->a[j]) >=  0.) continue;
376 #else
377         if (a->a[j] >=  0.) continue;
378 #endif
379         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
380       }
381     }
382     color = DRAW_CYAN;
383     for ( i=0; i<m; i++ ) {
384       y_l = m - i - 1.0; y_r = y_l + 1.0;
385       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
386         x_l = a->j[j] + shift; x_r = x_l + 1.0;
387         if (a->a[j] !=  0.) continue;
388         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
389       }
390     }
391     color = DRAW_RED;
392     for ( i=0; i<m; i++ ) {
393       y_l = m - i - 1.0; y_r = y_l + 1.0;
394       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
395         x_l = a->j[j] + shift; x_r = x_l + 1.0;
396 #if defined(PETSC_COMPLEX)
397         if (real(a->a[j]) <=  0.) continue;
398 #else
399         if (a->a[j] <=  0.) continue;
400 #endif
401         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
402       }
403     }
404     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
405   }
406   return 0;
407 }
408 
409 static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
410 {
411   Mat         A = (Mat) obj;
412   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
413   ViewerType  vtype;
414   int         ierr;
415 
416   if (!viewer) {
417     viewer = STDOUT_VIEWER_SELF;
418   }
419   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
420   if (vtype == MATLAB_VIEWER) {
421     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
422   }
423   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
424     return MatView_SeqAIJ_ASCII(A,viewer);
425   }
426   else if (vtype == BINARY_FILE_VIEWER) {
427     return MatView_SeqAIJ_Binary(A,viewer);
428   }
429   else if (vtype == DRAW_VIEWER) {
430     return MatView_SeqAIJ_Draw(A,viewer);
431   }
432   return 0;
433 }
434 
435 extern int Mat_AIJ_CheckInode(Mat);
436 static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
437 {
438   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
439   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
440   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
441   Scalar     *aa = a->a, *ap;
442 
443   if (mode == FLUSH_ASSEMBLY) return 0;
444 
445   for ( i=1; i<m; i++ ) {
446     /* move each row back by the amount of empty slots (fshift) before it*/
447     fshift += imax[i-1] - ailen[i-1];
448     if (fshift) {
449       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
450       N = ailen[i];
451       for ( j=0; j<N; j++ ) {
452         ip[j-fshift] = ip[j];
453         ap[j-fshift] = ap[j];
454       }
455     }
456     ai[i] = ai[i-1] + ailen[i-1];
457   }
458   if (m) {
459     fshift += imax[m-1] - ailen[m-1];
460     ai[m] = ai[m-1] + ailen[m-1];
461   }
462   /* reset ilen and imax for each row */
463   for ( i=0; i<m; i++ ) {
464     ailen[i] = imax[i] = ai[i+1] - ai[i];
465   }
466   a->nz = ai[m] + shift;
467 
468   /* diagonals may have moved, so kill the diagonal pointers */
469   if (fshift && a->diag) {
470     PetscFree(a->diag);
471     PLogObjectMemory(A,-(m+1)*sizeof(int));
472     a->diag = 0;
473   }
474   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n",
475            fshift,a->nz,m);
476   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n",
477            a->reallocs);
478   /* check out for identical nodes. If found, use inode functions */
479   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
480   return 0;
481 }
482 
483 static int MatZeroEntries_SeqAIJ(Mat A)
484 {
485   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
486   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
487   return 0;
488 }
489 
490 int MatDestroy_SeqAIJ(PetscObject obj)
491 {
492   Mat        A  = (Mat) obj;
493   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
494 
495 #if defined(PETSC_LOG)
496   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
497 #endif
498   PetscFree(a->a);
499   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
500   if (a->diag) PetscFree(a->diag);
501   if (a->ilen) PetscFree(a->ilen);
502   if (a->imax) PetscFree(a->imax);
503   if (a->solve_work) PetscFree(a->solve_work);
504   if (a->inode.size) PetscFree(a->inode.size);
505   PetscFree(a);
506   PLogObjectDestroy(A);
507   PetscHeaderDestroy(A);
508   return 0;
509 }
510 
511 static int MatCompress_SeqAIJ(Mat A)
512 {
513   return 0;
514 }
515 
516 static int MatSetOption_SeqAIJ(Mat A,MatOption op)
517 {
518   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
519   if      (op == ROW_ORIENTED)              a->roworiented = 1;
520   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
521   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
522   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
523   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
524   else if (op == ROWS_SORTED ||
525            op == SYMMETRIC_MATRIX ||
526            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
527            op == YES_NEW_DIAGONALS)
528     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
529   else if (op == NO_NEW_DIAGONALS)
530     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
531   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
532   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
533   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
534   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
535   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
536   else
537     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
538   return 0;
539 }
540 
541 static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
542 {
543   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
544   int        i,j, n,shift = a->indexshift;
545   Scalar     *x, zero = 0.0;
546 
547   VecSet(&zero,v);
548   VecGetArray(v,&x); VecGetLocalSize(v,&n);
549   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
550   for ( i=0; i<a->m; i++ ) {
551     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
552       if (a->j[j]+shift == i) {
553         x[i] = a->a[j];
554         break;
555       }
556     }
557   }
558   return 0;
559 }
560 
561 /* -------------------------------------------------------*/
562 /* Should check that shapes of vectors and matrices match */
563 /* -------------------------------------------------------*/
564 static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
565 {
566   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
567   Scalar     *x, *y, *v, alpha;
568   int        m = a->m, n, i, *idx, shift = a->indexshift;
569 
570   VecGetArray(xx,&x); VecGetArray(yy,&y);
571   PetscMemzero(y,a->n*sizeof(Scalar));
572   y = y + shift; /* shift for Fortran start by 1 indexing */
573   for ( i=0; i<m; i++ ) {
574     idx   = a->j + a->i[i] + shift;
575     v     = a->a + a->i[i] + shift;
576     n     = a->i[i+1] - a->i[i];
577     alpha = x[i];
578     while (n-->0) {y[*idx++] += alpha * *v++;}
579   }
580   PLogFlops(2*a->nz - a->n);
581   return 0;
582 }
583 
584 static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
585 {
586   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
587   Scalar     *x, *y, *v, alpha;
588   int        m = a->m, n, i, *idx,shift = a->indexshift;
589 
590   VecGetArray(xx,&x); VecGetArray(yy,&y);
591   if (zz != yy) VecCopy(zz,yy);
592   y = y + shift; /* shift for Fortran start by 1 indexing */
593   for ( i=0; i<m; i++ ) {
594     idx   = a->j + a->i[i] + shift;
595     v     = a->a + a->i[i] + shift;
596     n     = a->i[i+1] - a->i[i];
597     alpha = x[i];
598     while (n-->0) {y[*idx++] += alpha * *v++;}
599   }
600   return 0;
601 }
602 
603 static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
604 {
605   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
606   Scalar     *x, *y, *v, sum;
607   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
608 
609   VecGetArray(xx,&x); VecGetArray(yy,&y);
610   x    = x + shift; /* shift for Fortran start by 1 indexing */
611   idx  = a->j;
612   v    = a->a;
613   ii   = a->i;
614   for ( i=0; i<m; i++ ) {
615     n    = ii[1] - ii[0]; ii++;
616     sum  = 0.0;
617     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
618     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
619     while (n--) sum += *v++ * x[*idx++];
620     y[i] = sum;
621   }
622   PLogFlops(2*a->nz - m);
623   return 0;
624 }
625 
626 static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
627 {
628   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
629   Scalar     *x, *y, *z, *v, sum;
630   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
631 
632   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
633   x    = x + shift; /* shift for Fortran start by 1 indexing */
634   idx  = a->j;
635   v    = a->a;
636   ii   = a->i;
637   for ( i=0; i<m; i++ ) {
638     n    = ii[1] - ii[0]; ii++;
639     sum  = y[i];
640     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
641     while (n--) sum += *v++ * x[*idx++];
642     z[i] = sum;
643   }
644   PLogFlops(2*a->nz);
645   return 0;
646 }
647 
648 /*
649      Adds diagonal pointers to sparse matrix structure.
650 */
651 
652 int MatMarkDiag_SeqAIJ(Mat A)
653 {
654   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
655   int        i,j, *diag, m = a->m,shift = a->indexshift;
656 
657   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
658   PLogObjectMemory(A,(m+1)*sizeof(int));
659   for ( i=0; i<a->m; i++ ) {
660     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
661       if (a->j[j]+shift == i) {
662         diag[i] = j - shift;
663         break;
664       }
665     }
666   }
667   a->diag = diag;
668   return 0;
669 }
670 
671 static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
672                            double fshift,int its,Vec xx)
673 {
674   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
675   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
676   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
677 
678   VecGetArray(xx,&x); VecGetArray(bb,&b);
679   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
680   diag = a->diag;
681   xs   = x + shift; /* shifted by one for index start of a or a->j*/
682   if (flag == SOR_APPLY_UPPER) {
683    /* apply ( U + D/omega) to the vector */
684     bs = b + shift;
685     for ( i=0; i<m; i++ ) {
686         d    = fshift + a->a[diag[i] + shift];
687         n    = a->i[i+1] - diag[i] - 1;
688         idx  = a->j + diag[i] + (!shift);
689         v    = a->a + diag[i] + (!shift);
690         sum  = b[i]*d/omega;
691         SPARSEDENSEDOT(sum,bs,v,idx,n);
692         x[i] = sum;
693     }
694     return 0;
695   }
696   if (flag == SOR_APPLY_LOWER) {
697     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
698   }
699   else if (flag & SOR_EISENSTAT) {
700     /* Let  A = L + U + D; where L is lower trianglar,
701     U is upper triangular, E is diagonal; This routine applies
702 
703             (L + E)^{-1} A (U + E)^{-1}
704 
705     to a vector efficiently using Eisenstat's trick. This is for
706     the case of SSOR preconditioner, so E is D/omega where omega
707     is the relaxation factor.
708     */
709     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
710     scale = (2.0/omega) - 1.0;
711 
712     /*  x = (E + U)^{-1} b */
713     for ( i=m-1; i>=0; i-- ) {
714       d    = fshift + a->a[diag[i] + shift];
715       n    = a->i[i+1] - diag[i] - 1;
716       idx  = a->j + diag[i] + (!shift);
717       v    = a->a + diag[i] + (!shift);
718       sum  = b[i];
719       SPARSEDENSEMDOT(sum,xs,v,idx,n);
720       x[i] = omega*(sum/d);
721     }
722 
723     /*  t = b - (2*E - D)x */
724     v = a->a;
725     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
726 
727     /*  t = (E + L)^{-1}t */
728     ts = t + shift; /* shifted by one for index start of a or a->j*/
729     diag = a->diag;
730     for ( i=0; i<m; i++ ) {
731       d    = fshift + a->a[diag[i]+shift];
732       n    = diag[i] - a->i[i];
733       idx  = a->j + a->i[i] + shift;
734       v    = a->a + a->i[i] + shift;
735       sum  = t[i];
736       SPARSEDENSEMDOT(sum,ts,v,idx,n);
737       t[i] = omega*(sum/d);
738     }
739 
740     /*  x = x + t */
741     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
742     PetscFree(t);
743     return 0;
744   }
745   if (flag & SOR_ZERO_INITIAL_GUESS) {
746     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
747       for ( i=0; i<m; i++ ) {
748         d    = fshift + a->a[diag[i]+shift];
749         n    = diag[i] - a->i[i];
750         idx  = a->j + a->i[i] + shift;
751         v    = a->a + a->i[i] + shift;
752         sum  = b[i];
753         SPARSEDENSEMDOT(sum,xs,v,idx,n);
754         x[i] = omega*(sum/d);
755       }
756       xb = x;
757     }
758     else xb = b;
759     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
760         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
761       for ( i=0; i<m; i++ ) {
762         x[i] *= a->a[diag[i]+shift];
763       }
764     }
765     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
766       for ( i=m-1; i>=0; i-- ) {
767         d    = fshift + a->a[diag[i] + shift];
768         n    = a->i[i+1] - diag[i] - 1;
769         idx  = a->j + diag[i] + (!shift);
770         v    = a->a + diag[i] + (!shift);
771         sum  = xb[i];
772         SPARSEDENSEMDOT(sum,xs,v,idx,n);
773         x[i] = omega*(sum/d);
774       }
775     }
776     its--;
777   }
778   while (its--) {
779     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
780       for ( i=0; i<m; i++ ) {
781         d    = fshift + a->a[diag[i]+shift];
782         n    = a->i[i+1] - a->i[i];
783         idx  = a->j + a->i[i] + shift;
784         v    = a->a + a->i[i] + shift;
785         sum  = b[i];
786         SPARSEDENSEMDOT(sum,xs,v,idx,n);
787         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
788       }
789     }
790     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
791       for ( i=m-1; i>=0; i-- ) {
792         d    = fshift + a->a[diag[i] + shift];
793         n    = a->i[i+1] - a->i[i];
794         idx  = a->j + a->i[i] + shift;
795         v    = a->a + a->i[i] + shift;
796         sum  = b[i];
797         SPARSEDENSEMDOT(sum,xs,v,idx,n);
798         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
799       }
800     }
801   }
802   return 0;
803 }
804 
805 static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
806 {
807   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
808   if (nz)      *nz      = a->nz;
809   if (nzalloc) *nzalloc = a->maxnz;
810   if (mem)     *mem     = (int)A->mem;
811   return 0;
812 }
813 
814 extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
815 extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
816 extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
817 extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
818 extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
819 extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
820 extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
821 
822 static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
823 {
824   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
825   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
826 
827   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
828   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
829   if (diag) {
830     for ( i=0; i<N; i++ ) {
831       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
832       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
833         a->ilen[rows[i]] = 1;
834         a->a[a->i[rows[i]]+shift] = *diag;
835         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
836       }
837       else {
838         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
839         CHKERRQ(ierr);
840       }
841     }
842   }
843   else {
844     for ( i=0; i<N; i++ ) {
845       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
846       a->ilen[rows[i]] = 0;
847     }
848   }
849   ISRestoreIndices(is,&rows);
850   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
851   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
852   return 0;
853 }
854 
855 static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
856 {
857   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
858   *m = a->m; *n = a->n;
859   return 0;
860 }
861 
862 static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
863 {
864   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
865   *m = 0; *n = a->m;
866   return 0;
867 }
868 int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
869 {
870   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
871   int        *itmp,i,shift = a->indexshift;
872 
873   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
874 
875   *nz = a->i[row+1] - a->i[row];
876   if (v) *v = a->a + a->i[row] + shift;
877   if (idx) {
878     itmp = a->j + a->i[row] + shift;
879     if (*nz && shift) {
880       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
881       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
882     } else if (*nz) {
883       *idx = itmp;
884     }
885     else *idx = 0;
886   }
887   return 0;
888 }
889 
890 int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
891 {
892   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
893   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
894   return 0;
895 }
896 
897 static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
898 {
899   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
900   Scalar     *v = a->a;
901   double     sum = 0.0;
902   int        i, j,shift = a->indexshift;
903 
904   if (type == NORM_FROBENIUS) {
905     for (i=0; i<a->nz; i++ ) {
906 #if defined(PETSC_COMPLEX)
907       sum += real(conj(*v)*(*v)); v++;
908 #else
909       sum += (*v)*(*v); v++;
910 #endif
911     }
912     *norm = sqrt(sum);
913   }
914   else if (type == NORM_1) {
915     double *tmp;
916     int    *jj = a->j;
917     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
918     PetscMemzero(tmp,a->n*sizeof(double));
919     *norm = 0.0;
920     for ( j=0; j<a->nz; j++ ) {
921         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
922     }
923     for ( j=0; j<a->n; j++ ) {
924       if (tmp[j] > *norm) *norm = tmp[j];
925     }
926     PetscFree(tmp);
927   }
928   else if (type == NORM_INFINITY) {
929     *norm = 0.0;
930     for ( j=0; j<a->m; j++ ) {
931       v = a->a + a->i[j] + shift;
932       sum = 0.0;
933       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
934         sum += PetscAbsScalar(*v); v++;
935       }
936       if (sum > *norm) *norm = sum;
937     }
938   }
939   else {
940     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
941   }
942   return 0;
943 }
944 
945 static int MatTranspose_SeqAIJ(Mat A,Mat *B)
946 {
947   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
948   Mat        C;
949   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
950   int        shift = a->indexshift;
951   Scalar     *array = a->a;
952 
953   if (B == PETSC_NULL && m != a->n)
954     SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place");
955   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
956   PetscMemzero(col,(1+a->n)*sizeof(int));
957   if (shift) {
958     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
959   }
960   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
961   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
962   PetscFree(col);
963   for ( i=0; i<m; i++ ) {
964     len = ai[i+1]-ai[i];
965     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
966     array += len; aj += len;
967   }
968   if (shift) {
969     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
970   }
971 
972   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
973   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
974 
975   if (B != PETSC_NULL) {
976     *B = C;
977   } else {
978     /* This isn't really an in-place transpose */
979     PetscFree(a->a);
980     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
981     if (a->diag) PetscFree(a->diag);
982     if (a->ilen) PetscFree(a->ilen);
983     if (a->imax) PetscFree(a->imax);
984     if (a->solve_work) PetscFree(a->solve_work);
985     if (a->inode.size) PetscFree(a->inode.size);
986     PetscFree(a);
987     PetscMemcpy(A,C,sizeof(struct _Mat));
988     PetscHeaderDestroy(C);
989   }
990   return 0;
991 }
992 
993 static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
994 {
995   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
996   Scalar     *l,*r,x,*v;
997   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
998 
999   if (ll) {
1000     VecGetArray(ll,&l); VecGetSize(ll,&m);
1001     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
1002     v = a->a;
1003     for ( i=0; i<m; i++ ) {
1004       x = l[i];
1005       M = a->i[i+1] - a->i[i];
1006       for ( j=0; j<M; j++ ) { (*v++) *= x;}
1007     }
1008   }
1009   if (rr) {
1010     VecGetArray(rr,&r); VecGetSize(rr,&n);
1011     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
1012     v = a->a; jj = a->j;
1013     for ( i=0; i<nz; i++ ) {
1014       (*v++) *= r[*jj++ + shift];
1015     }
1016   }
1017   return 0;
1018 }
1019 
1020 static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
1021 {
1022   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1023   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
1024   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1025   register int sum,lensi;
1026   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
1027   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
1028   Scalar       *a_new,*mat_a;
1029   Mat          C;
1030 
1031   ierr = ISSorted(iscol,(PetscTruth*)&i);
1032   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IS is not sorted");
1033 
1034   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
1035   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
1036   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
1037 
1038   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
1039     /* special case of contiguous rows */
1040     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
1041     starts = lens + ncols;
1042     /* loop over new rows determining lens and starting points */
1043     for (i=0; i<nrows; i++) {
1044       kstart  = ai[irow[i]]+shift;
1045       kend    = kstart + ailen[irow[i]];
1046       for ( k=kstart; k<kend; k++ ) {
1047         if (aj[k]+shift >= first) {
1048           starts[i] = k;
1049           break;
1050 	}
1051       }
1052       sum = 0;
1053       while (k < kend) {
1054         if (aj[k++]+shift >= first+ncols) break;
1055         sum++;
1056       }
1057       lens[i] = sum;
1058     }
1059     /* create submatrix */
1060     if (scall == MAT_REUSE_MATRIX) {
1061       int n_cols,n_rows;
1062       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1063       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1064       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
1065       C = *B;
1066     }
1067     else {
1068       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
1069     }
1070     c = (Mat_SeqAIJ*) C->data;
1071 
1072     /* loop over rows inserting into submatrix */
1073     a_new    = c->a;
1074     j_new    = c->j;
1075     i_new    = c->i;
1076     i_new[0] = -shift;
1077     for (i=0; i<nrows; i++) {
1078       ii    = starts[i];
1079       lensi = lens[i];
1080       for ( k=0; k<lensi; k++ ) {
1081         *j_new++ = aj[ii+k] - first;
1082       }
1083       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1084       a_new      += lensi;
1085       i_new[i+1]  = i_new[i] + lensi;
1086       c->ilen[i]  = lensi;
1087     }
1088     PetscFree(lens);
1089   }
1090   else {
1091     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
1092     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
1093     ssmap = smap + shift;
1094     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1095     PetscMemzero(smap,oldcols*sizeof(int));
1096     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
1097     /* determine lens of each row */
1098     for (i=0; i<nrows; i++) {
1099       kstart  = ai[irow[i]]+shift;
1100       kend    = kstart + a->ilen[irow[i]];
1101       lens[i] = 0;
1102       for ( k=kstart; k<kend; k++ ) {
1103         if (ssmap[aj[k]]) {
1104           lens[i]++;
1105         }
1106       }
1107     }
1108     /* Create and fill new matrix */
1109     if (scall == MAT_REUSE_MATRIX) {
1110       c = (Mat_SeqAIJ *)((*B)->data);
1111 
1112       if (c->m  != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1113       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1114         SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros");
1115       }
1116       PetscMemzero(c->ilen,c->m*sizeof(int));
1117       C = *B;
1118     }
1119     else {
1120       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
1121     }
1122     c = (Mat_SeqAIJ *)(C->data);
1123     for (i=0; i<nrows; i++) {
1124       row    = irow[i];
1125       nznew  = 0;
1126       kstart = ai[row]+shift;
1127       kend   = kstart + a->ilen[row];
1128       mat_i  = c->i[i]+shift;
1129       mat_j  = c->j + mat_i;
1130       mat_a  = c->a + mat_i;
1131       mat_ilen = c->ilen + i;
1132       for ( k=kstart; k<kend; k++ ) {
1133         if ((tcol=ssmap[a->j[k]])) {
1134           *mat_j++ = tcol - (!shift);
1135           *mat_a++ = a->a[k];
1136           (*mat_ilen)++;
1137 
1138         }
1139       }
1140     }
1141     /* Free work space */
1142     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1143     PetscFree(smap); PetscFree(lens);
1144   }
1145   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1146   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1147 
1148   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1149   *B = C;
1150   return 0;
1151 }
1152 
1153 /*
1154      note: This can only work for identity for row and col. It would
1155    be good to check this and otherwise generate an error.
1156 */
1157 static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1158 {
1159   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1160   int        ierr;
1161   Mat        outA;
1162 
1163   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1164 
1165   outA          = inA;
1166   inA->factor   = FACTOR_LU;
1167   a->row        = row;
1168   a->col        = col;
1169 
1170   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
1171 
1172   if (!a->diag) {
1173     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
1174   }
1175   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1176   return 0;
1177 }
1178 
1179 #include "pinclude/plapack.h"
1180 static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1181 {
1182   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1183   int        one = 1;
1184   BLscal_( &a->nz, alpha, a->a, &one );
1185   PLogFlops(a->nz);
1186   return 0;
1187 }
1188 
1189 static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1190                                     Mat **B)
1191 {
1192   int ierr,i;
1193 
1194   if (scall == MAT_INITIAL_MATRIX) {
1195     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1196   }
1197 
1198   for ( i=0; i<n; i++ ) {
1199     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1200   }
1201   return 0;
1202 }
1203 
1204 static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
1205 {
1206   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1207   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
1208   int        start, end, *ai, *aj;
1209   char       *table;
1210   shift = a->indexshift;
1211   m     = a->m;
1212   ai    = a->i;
1213   aj    = a->j+shift;
1214 
1215   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
1216 
1217   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
1218   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1219 
1220   for ( i=0; i<is_max; i++ ) {
1221     /* Initialise the two local arrays */
1222     isz  = 0;
1223     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1224 
1225                 /* Extract the indices, assume there can be duplicate entries */
1226     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
1227     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1228 
1229     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1230     for ( j=0; j<n ; ++j){
1231       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
1232     }
1233     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
1234     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1235 
1236     k = 0;
1237     for ( j=0; j<ov; j++){ /* for each overlap*/
1238       n = isz;
1239       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1240         row   = nidx[k];
1241         start = ai[row];
1242         end   = ai[row+1];
1243         for ( l = start; l<end ; l++){
1244           val = aj[l] + shift;
1245           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1246         }
1247       }
1248     }
1249     ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1250   }
1251   PetscFree(table);
1252   PetscFree(nidx);
1253   return 0;
1254 }
1255 
1256 int MatPrintHelp_SeqAIJ(Mat A)
1257 {
1258   static int called = 0;
1259   MPI_Comm   comm = A->comm;
1260 
1261   if (called) return 0; else called = 1;
1262   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1263   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
1264   PetscPrintf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
1265   PetscPrintf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
1266   PetscPrintf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1267 #if defined(HAVE_ESSL)
1268   PetscPrintf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1269 #endif
1270   return 0;
1271 }
1272 static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1273 /* -------------------------------------------------------------------*/
1274 static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
1275        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1276        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1277        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
1278        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
1279        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
1280        MatLUFactor_SeqAIJ,0,
1281        MatRelax_SeqAIJ,
1282        MatTranspose_SeqAIJ,
1283        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1284        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
1285        0,MatAssemblyEnd_SeqAIJ,
1286        MatCompress_SeqAIJ,
1287        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
1288        MatGetReordering_SeqAIJ,
1289        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
1290        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
1291        MatILUFactorSymbolic_SeqAIJ,0,
1292        0,0,MatConvert_SeqAIJ,
1293        MatGetSubMatrix_SeqAIJ,0,
1294        MatConvertSameType_SeqAIJ,0,0,
1295        MatILUFactor_SeqAIJ,0,0,
1296        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1297        MatGetValues_SeqAIJ,0,
1298        MatPrintHelp_SeqAIJ,
1299        MatScale_SeqAIJ};
1300 
1301 extern int MatUseSuperLU_SeqAIJ(Mat);
1302 extern int MatUseEssl_SeqAIJ(Mat);
1303 extern int MatUseDXML_SeqAIJ(Mat);
1304 
1305 /*@C
1306    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
1307    (the default parallel PETSc format).  For good matrix assembly performance
1308    the user should preallocate the matrix storage by setting the parameter nz
1309    (or nzz).  By setting these parameters accurately, performance can be
1310    increased by more than a factor of 50.
1311 
1312    Input Parameters:
1313 .  comm - MPI communicator, set to MPI_COMM_SELF
1314 .  m - number of rows
1315 .  n - number of columns
1316 .  nz - number of nonzeros per row (same for all rows)
1317 .  nzz - number of nonzeros per row or null (possibly different for each row)
1318 
1319    Output Parameter:
1320 .  A - the matrix
1321 
1322    Notes:
1323    The AIJ format (also called the Yale sparse matrix format or
1324    compressed row storage), is fully compatible with standard Fortran 77
1325    storage.  That is, the stored row and column indices can begin at
1326    either one (as in Fortran) or zero.  See the users manual for details.
1327 
1328    Specify the preallocated storage with either nz or nnz (not both).
1329    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1330    allocation.  For additional details, see the users manual chapter on
1331    matrices and the file $(PETSC_DIR)/Performance.
1332 
1333    By default, this format uses inodes (identical nodes) when possible, to
1334    improve numerical efficiency of Matrix vector products and solves. We
1335    search for consecutive rows with the same nonzero structure, thereby
1336    reusing matrix information to achieve increased efficiency.
1337 
1338    Options Database Keys:
1339 $    -mat_aij_no_inode  - Do not use inodes
1340 $    -mat_aij_inode_limit <limit> - Set inode limit.
1341 $        (max limit=5)
1342 $    -mat_aij_oneindex - Internally use indexing starting at 1
1343 $        rather than 0.  Note: When calling MatSetValues(),
1344 $        the user still MUST index entries starting at 0!
1345 
1346 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
1347 @*/
1348 int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
1349 {
1350   Mat        B;
1351   Mat_SeqAIJ *b;
1352   int        i,len,ierr, flg;
1353 
1354   *A      = 0;
1355   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1356   PLogObjectCreate(B);
1357   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1358   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1359   B->destroy          = MatDestroy_SeqAIJ;
1360   B->view             = MatView_SeqAIJ;
1361   B->factor           = 0;
1362   B->lupivotthreshold = 1.0;
1363   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
1364                           &flg); CHKERRQ(ierr);
1365   b->ilu_preserve_row_sums = PETSC_FALSE;
1366   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
1367                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1368   b->row              = 0;
1369   b->col              = 0;
1370   b->indexshift       = 0;
1371   b->reallocs         = 0;
1372   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
1373   if (flg) b->indexshift = -1;
1374 
1375   b->m       = m;
1376   b->n       = n;
1377   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1378   if (nnz == PETSC_NULL) {
1379     if (nz == PETSC_DEFAULT) nz = 10;
1380     else if (nz <= 0)        nz = 1;
1381     for ( i=0; i<m; i++ ) b->imax[i] = nz;
1382     nz = nz*m;
1383   }
1384   else {
1385     nz = 0;
1386     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
1387   }
1388 
1389   /* allocate the matrix space */
1390   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
1391   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1392   b->j  = (int *) (b->a + nz);
1393   PetscMemzero(b->j,nz*sizeof(int));
1394   b->i  = b->j + nz;
1395   b->singlemalloc = 1;
1396 
1397   b->i[0] = -b->indexshift;
1398   for (i=1; i<m+1; i++) {
1399     b->i[i] = b->i[i-1] + b->imax[i-1];
1400   }
1401 
1402   /* b->ilen will count nonzeros in each row so far. */
1403   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1404   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1405   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
1406 
1407   b->nz               = 0;
1408   b->maxnz            = nz;
1409   b->sorted           = 0;
1410   b->roworiented      = 1;
1411   b->nonew            = 0;
1412   b->diag             = 0;
1413   b->solve_work       = 0;
1414   b->spptr            = 0;
1415   b->inode.node_count = 0;
1416   b->inode.size       = 0;
1417   b->inode.limit      = 5;
1418   b->inode.max_limit  = 5;
1419 
1420   *A = B;
1421   /*  SuperLU is not currently supported through PETSc */
1422 #if defined(HAVE_SUPERLU)
1423   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
1424   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
1425 #endif
1426   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
1427   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
1428   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
1429   if (flg) {
1430     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1431     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
1432   }
1433   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1434   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1435   return 0;
1436 }
1437 
1438 int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
1439 {
1440   Mat        C;
1441   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1442   int        i,len, m = a->m,shift = a->indexshift;
1443 
1444   *B = 0;
1445   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1446   PLogObjectCreate(C);
1447   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1448   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1449   C->destroy    = MatDestroy_SeqAIJ;
1450   C->view       = MatView_SeqAIJ;
1451   C->factor     = A->factor;
1452   c->row        = 0;
1453   c->col        = 0;
1454   c->indexshift = shift;
1455   C->assembled  = PETSC_TRUE;
1456 
1457   c->m          = a->m;
1458   c->n          = a->n;
1459 
1460   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1461   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
1462   for ( i=0; i<m; i++ ) {
1463     c->imax[i] = a->imax[i];
1464     c->ilen[i] = a->ilen[i];
1465   }
1466 
1467   /* allocate the matrix space */
1468   c->singlemalloc = 1;
1469   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1470   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1471   c->j  = (int *) (c->a + a->i[m] + shift);
1472   c->i  = c->j + a->i[m] + shift;
1473   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
1474   if (m > 0) {
1475     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1476     if (cpvalues == COPY_VALUES) {
1477       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
1478     }
1479   }
1480 
1481   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1482   c->sorted      = a->sorted;
1483   c->roworiented = a->roworiented;
1484   c->nonew       = a->nonew;
1485   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
1486 
1487   if (a->diag) {
1488     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1489     PLogObjectMemory(C,(m+1)*sizeof(int));
1490     for ( i=0; i<m; i++ ) {
1491       c->diag[i] = a->diag[i];
1492     }
1493   }
1494   else c->diag          = 0;
1495   c->inode.limit        = a->inode.limit;
1496   c->inode.max_limit    = a->inode.max_limit;
1497   if (a->inode.size){
1498     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1499     c->inode.node_count = a->inode.node_count;
1500     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1501   } else {
1502     c->inode.size       = 0;
1503     c->inode.node_count = 0;
1504   }
1505   c->nz                 = a->nz;
1506   c->maxnz              = a->maxnz;
1507   c->solve_work         = 0;
1508   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1509 
1510   *B = C;
1511   return 0;
1512 }
1513 
1514 int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
1515 {
1516   Mat_SeqAIJ   *a;
1517   Mat          B;
1518   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1519   MPI_Comm     comm;
1520 
1521   PetscObjectGetComm((PetscObject) viewer,&comm);
1522   MPI_Comm_size(comm,&size);
1523   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
1524   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
1525   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1526   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
1527   M = header[1]; N = header[2]; nz = header[3];
1528 
1529   /* read in row lengths */
1530   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1531   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
1532 
1533   /* create our matrix */
1534   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1535   B = *A;
1536   a = (Mat_SeqAIJ *) B->data;
1537   shift = a->indexshift;
1538 
1539   /* read in column indices and adjust for Fortran indexing*/
1540   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
1541   if (shift) {
1542     for ( i=0; i<nz; i++ ) {
1543       a->j[i] += 1;
1544     }
1545   }
1546 
1547   /* read in nonzero values */
1548   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
1549 
1550   /* set matrix "i" values */
1551   a->i[0] = -shift;
1552   for ( i=1; i<= M; i++ ) {
1553     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1554     a->ilen[i-1] = rowlengths[i-1];
1555   }
1556   PetscFree(rowlengths);
1557 
1558   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1559   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1560   return 0;
1561 }
1562 
1563 static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
1564 {
1565   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
1566 
1567   if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type");
1568 
1569   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
1570   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1571       (a->indexshift != b->indexshift)) {
1572     *flg = PETSC_FALSE; return 0;
1573   }
1574 
1575   /* if the a->i are the same */
1576   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
1577     *flg = PETSC_FALSE; return 0;
1578   }
1579 
1580   /* if a->j are the same */
1581   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
1582     *flg = PETSC_FALSE; return 0;
1583   }
1584 
1585   /* if a->a are the same */
1586   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
1587     *flg = PETSC_FALSE; return 0;
1588   }
1589   *flg = PETSC_TRUE;
1590   return 0;
1591 
1592 }
1593