xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 17f4b42bc68283ebfe06c2f879019e653f788d8f)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: aij.c,v 1.163 1996/03/26 01:01:27 bsmith Exp bsmith $";
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   return 0;
507 }
508 
509 static int MatCompress_SeqAIJ(Mat A)
510 {
511   return 0;
512 }
513 
514 static int MatSetOption_SeqAIJ(Mat A,MatOption op)
515 {
516   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
517   if      (op == ROW_ORIENTED)              a->roworiented = 1;
518   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
519   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
520   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
521   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
522   else if (op == ROWS_SORTED ||
523            op == SYMMETRIC_MATRIX ||
524            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
525            op == YES_NEW_DIAGONALS)
526     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
527   else if (op == NO_NEW_DIAGONALS)
528     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
529   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
530   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
531   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
532   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
533   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
534   else
535     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
536   return 0;
537 }
538 
539 static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
540 {
541   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
542   int        i,j, n,shift = a->indexshift;
543   Scalar     *x, zero = 0.0;
544 
545   VecSet(&zero,v);
546   VecGetArray(v,&x); VecGetLocalSize(v,&n);
547   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
548   for ( i=0; i<a->m; i++ ) {
549     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
550       if (a->j[j]+shift == i) {
551         x[i] = a->a[j];
552         break;
553       }
554     }
555   }
556   return 0;
557 }
558 
559 /* -------------------------------------------------------*/
560 /* Should check that shapes of vectors and matrices match */
561 /* -------------------------------------------------------*/
562 static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
563 {
564   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
565   Scalar     *x, *y, *v, alpha;
566   int        m = a->m, n, i, *idx, shift = a->indexshift;
567 
568   VecGetArray(xx,&x); VecGetArray(yy,&y);
569   PetscMemzero(y,a->n*sizeof(Scalar));
570   y = y + shift; /* shift for Fortran start by 1 indexing */
571   for ( i=0; i<m; i++ ) {
572     idx   = a->j + a->i[i] + shift;
573     v     = a->a + a->i[i] + shift;
574     n     = a->i[i+1] - a->i[i];
575     alpha = x[i];
576     while (n-->0) {y[*idx++] += alpha * *v++;}
577   }
578   PLogFlops(2*a->nz - a->n);
579   return 0;
580 }
581 
582 static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
583 {
584   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
585   Scalar     *x, *y, *v, alpha;
586   int        m = a->m, n, i, *idx,shift = a->indexshift;
587 
588   VecGetArray(xx,&x); VecGetArray(yy,&y);
589   if (zz != yy) VecCopy(zz,yy);
590   y = y + shift; /* shift for Fortran start by 1 indexing */
591   for ( i=0; i<m; i++ ) {
592     idx   = a->j + a->i[i] + shift;
593     v     = a->a + a->i[i] + shift;
594     n     = a->i[i+1] - a->i[i];
595     alpha = x[i];
596     while (n-->0) {y[*idx++] += alpha * *v++;}
597   }
598   return 0;
599 }
600 
601 static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
602 {
603   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
604   Scalar     *x, *y, *v, sum;
605   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
606 
607   VecGetArray(xx,&x); VecGetArray(yy,&y);
608   x    = x + shift; /* shift for Fortran start by 1 indexing */
609   idx  = a->j;
610   v    = a->a;
611   ii   = a->i;
612   for ( i=0; i<m; i++ ) {
613     n    = ii[1] - ii[0]; ii++;
614     sum  = 0.0;
615     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
616     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
617     while (n--) sum += *v++ * x[*idx++];
618     y[i] = sum;
619   }
620   PLogFlops(2*a->nz - m);
621   return 0;
622 }
623 
624 static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
625 {
626   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
627   Scalar     *x, *y, *z, *v, sum;
628   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
629 
630   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
631   x    = x + shift; /* shift for Fortran start by 1 indexing */
632   idx  = a->j;
633   v    = a->a;
634   ii   = a->i;
635   for ( i=0; i<m; i++ ) {
636     n    = ii[1] - ii[0]; ii++;
637     sum  = y[i];
638     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
639     while (n--) sum += *v++ * x[*idx++];
640     z[i] = sum;
641   }
642   PLogFlops(2*a->nz);
643   return 0;
644 }
645 
646 /*
647      Adds diagonal pointers to sparse matrix structure.
648 */
649 
650 int MatMarkDiag_SeqAIJ(Mat A)
651 {
652   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
653   int        i,j, *diag, m = a->m,shift = a->indexshift;
654 
655   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
656   PLogObjectMemory(A,(m+1)*sizeof(int));
657   for ( i=0; i<a->m; i++ ) {
658     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
659       if (a->j[j]+shift == i) {
660         diag[i] = j - shift;
661         break;
662       }
663     }
664   }
665   a->diag = diag;
666   return 0;
667 }
668 
669 static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
670                            double fshift,int its,Vec xx)
671 {
672   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
673   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
674   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
675 
676   VecGetArray(xx,&x); VecGetArray(bb,&b);
677   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
678   diag = a->diag;
679   xs   = x + shift; /* shifted by one for index start of a or a->j*/
680   if (flag == SOR_APPLY_UPPER) {
681    /* apply ( U + D/omega) to the vector */
682     bs = b + shift;
683     for ( i=0; i<m; i++ ) {
684         d    = fshift + a->a[diag[i] + shift];
685         n    = a->i[i+1] - diag[i] - 1;
686         idx  = a->j + diag[i] + (!shift);
687         v    = a->a + diag[i] + (!shift);
688         sum  = b[i]*d/omega;
689         SPARSEDENSEDOT(sum,bs,v,idx,n);
690         x[i] = sum;
691     }
692     return 0;
693   }
694   if (flag == SOR_APPLY_LOWER) {
695     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
696   }
697   else if (flag & SOR_EISENSTAT) {
698     /* Let  A = L + U + D; where L is lower trianglar,
699     U is upper triangular, E is diagonal; This routine applies
700 
701             (L + E)^{-1} A (U + E)^{-1}
702 
703     to a vector efficiently using Eisenstat's trick. This is for
704     the case of SSOR preconditioner, so E is D/omega where omega
705     is the relaxation factor.
706     */
707     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
708     scale = (2.0/omega) - 1.0;
709 
710     /*  x = (E + U)^{-1} b */
711     for ( i=m-1; i>=0; i-- ) {
712       d    = fshift + a->a[diag[i] + shift];
713       n    = a->i[i+1] - diag[i] - 1;
714       idx  = a->j + diag[i] + (!shift);
715       v    = a->a + diag[i] + (!shift);
716       sum  = b[i];
717       SPARSEDENSEMDOT(sum,xs,v,idx,n);
718       x[i] = omega*(sum/d);
719     }
720 
721     /*  t = b - (2*E - D)x */
722     v = a->a;
723     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
724 
725     /*  t = (E + L)^{-1}t */
726     ts = t + shift; /* shifted by one for index start of a or a->j*/
727     diag = a->diag;
728     for ( i=0; i<m; i++ ) {
729       d    = fshift + a->a[diag[i]+shift];
730       n    = diag[i] - a->i[i];
731       idx  = a->j + a->i[i] + shift;
732       v    = a->a + a->i[i] + shift;
733       sum  = t[i];
734       SPARSEDENSEMDOT(sum,ts,v,idx,n);
735       t[i] = omega*(sum/d);
736     }
737 
738     /*  x = x + t */
739     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
740     PetscFree(t);
741     return 0;
742   }
743   if (flag & SOR_ZERO_INITIAL_GUESS) {
744     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
745       for ( i=0; i<m; i++ ) {
746         d    = fshift + a->a[diag[i]+shift];
747         n    = diag[i] - a->i[i];
748         idx  = a->j + a->i[i] + shift;
749         v    = a->a + a->i[i] + shift;
750         sum  = b[i];
751         SPARSEDENSEMDOT(sum,xs,v,idx,n);
752         x[i] = omega*(sum/d);
753       }
754       xb = x;
755     }
756     else xb = b;
757     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
758         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
759       for ( i=0; i<m; i++ ) {
760         x[i] *= a->a[diag[i]+shift];
761       }
762     }
763     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
764       for ( i=m-1; i>=0; i-- ) {
765         d    = fshift + a->a[diag[i] + shift];
766         n    = a->i[i+1] - diag[i] - 1;
767         idx  = a->j + diag[i] + (!shift);
768         v    = a->a + diag[i] + (!shift);
769         sum  = xb[i];
770         SPARSEDENSEMDOT(sum,xs,v,idx,n);
771         x[i] = omega*(sum/d);
772       }
773     }
774     its--;
775   }
776   while (its--) {
777     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
778       for ( i=0; i<m; i++ ) {
779         d    = fshift + a->a[diag[i]+shift];
780         n    = a->i[i+1] - a->i[i];
781         idx  = a->j + a->i[i] + shift;
782         v    = a->a + a->i[i] + shift;
783         sum  = b[i];
784         SPARSEDENSEMDOT(sum,xs,v,idx,n);
785         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
786       }
787     }
788     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
789       for ( i=m-1; i>=0; i-- ) {
790         d    = fshift + a->a[diag[i] + shift];
791         n    = a->i[i+1] - a->i[i];
792         idx  = a->j + a->i[i] + shift;
793         v    = a->a + a->i[i] + shift;
794         sum  = b[i];
795         SPARSEDENSEMDOT(sum,xs,v,idx,n);
796         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
797       }
798     }
799   }
800   return 0;
801 }
802 
803 static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
804 {
805   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
806   if (nz)      *nz      = a->nz;
807   if (nzalloc) *nzalloc = a->maxnz;
808   if (mem)     *mem     = (int)A->mem;
809   return 0;
810 }
811 
812 extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
813 extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
814 extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
815 extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
816 extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
817 extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
818 extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
819 
820 static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
821 {
822   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
823   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
824 
825   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
826   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
827   if (diag) {
828     for ( i=0; i<N; i++ ) {
829       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
830       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
831         a->ilen[rows[i]] = 1;
832         a->a[a->i[rows[i]]+shift] = *diag;
833         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
834       }
835       else {
836         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
837         CHKERRQ(ierr);
838       }
839     }
840   }
841   else {
842     for ( i=0; i<N; i++ ) {
843       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
844       a->ilen[rows[i]] = 0;
845     }
846   }
847   ISRestoreIndices(is,&rows);
848   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
849   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
850   return 0;
851 }
852 
853 static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
854 {
855   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
856   *m = a->m; *n = a->n;
857   return 0;
858 }
859 
860 static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
861 {
862   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
863   *m = 0; *n = a->m;
864   return 0;
865 }
866 int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
867 {
868   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
869   int        *itmp,i,shift = a->indexshift;
870 
871   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
872 
873   *nz = a->i[row+1] - a->i[row];
874   if (v) *v = a->a + a->i[row] + shift;
875   if (idx) {
876     itmp = a->j + a->i[row] + shift;
877     if (*nz && shift) {
878       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
879       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
880     } else if (*nz) {
881       *idx = itmp;
882     }
883     else *idx = 0;
884   }
885   return 0;
886 }
887 
888 int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
889 {
890   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
891   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
892   return 0;
893 }
894 
895 static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
896 {
897   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
898   Scalar     *v = a->a;
899   double     sum = 0.0;
900   int        i, j,shift = a->indexshift;
901 
902   if (type == NORM_FROBENIUS) {
903     for (i=0; i<a->nz; i++ ) {
904 #if defined(PETSC_COMPLEX)
905       sum += real(conj(*v)*(*v)); v++;
906 #else
907       sum += (*v)*(*v); v++;
908 #endif
909     }
910     *norm = sqrt(sum);
911   }
912   else if (type == NORM_1) {
913     double *tmp;
914     int    *jj = a->j;
915     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
916     PetscMemzero(tmp,a->n*sizeof(double));
917     *norm = 0.0;
918     for ( j=0; j<a->nz; j++ ) {
919         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
920     }
921     for ( j=0; j<a->n; j++ ) {
922       if (tmp[j] > *norm) *norm = tmp[j];
923     }
924     PetscFree(tmp);
925   }
926   else if (type == NORM_INFINITY) {
927     *norm = 0.0;
928     for ( j=0; j<a->m; j++ ) {
929       v = a->a + a->i[j] + shift;
930       sum = 0.0;
931       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
932         sum += PetscAbsScalar(*v); v++;
933       }
934       if (sum > *norm) *norm = sum;
935     }
936   }
937   else {
938     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
939   }
940   return 0;
941 }
942 
943 static int MatTranspose_SeqAIJ(Mat A,Mat *B)
944 {
945   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
946   Mat        C;
947   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
948   int        shift = a->indexshift;
949   Scalar     *array = a->a;
950 
951   if (B == PETSC_NULL && m != a->n)
952     SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place");
953   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
954   PetscMemzero(col,(1+a->n)*sizeof(int));
955   if (shift) {
956     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
957   }
958   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
959   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
960   PetscFree(col);
961   for ( i=0; i<m; i++ ) {
962     len = ai[i+1]-ai[i];
963     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
964     array += len; aj += len;
965   }
966   if (shift) {
967     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
968   }
969 
970   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
971   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
972 
973   if (B != PETSC_NULL) {
974     *B = C;
975   } else {
976     /* This isn't really an in-place transpose */
977     PetscFree(a->a);
978     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
979     if (a->diag) PetscFree(a->diag);
980     if (a->ilen) PetscFree(a->ilen);
981     if (a->imax) PetscFree(a->imax);
982     if (a->solve_work) PetscFree(a->solve_work);
983     PetscFree(a);
984     PetscMemcpy(A,C,sizeof(struct _Mat));
985     PetscHeaderDestroy(C);
986   }
987   return 0;
988 }
989 
990 static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
991 {
992   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
993   Scalar     *l,*r,x,*v;
994   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
995 
996   if (ll) {
997     VecGetArray(ll,&l); VecGetSize(ll,&m);
998     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
999     v = a->a;
1000     for ( i=0; i<m; i++ ) {
1001       x = l[i];
1002       M = a->i[i+1] - a->i[i];
1003       for ( j=0; j<M; j++ ) { (*v++) *= x;}
1004     }
1005   }
1006   if (rr) {
1007     VecGetArray(rr,&r); VecGetSize(rr,&n);
1008     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
1009     v = a->a; jj = a->j;
1010     for ( i=0; i<nz; i++ ) {
1011       (*v++) *= r[*jj++ + shift];
1012     }
1013   }
1014   return 0;
1015 }
1016 
1017 static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
1018 {
1019   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1020   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
1021   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1022   register int sum,lensi;
1023   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
1024   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
1025   Scalar       *a_new,*mat_a;
1026   Mat          C;
1027 
1028   ierr = ISSorted(iscol,(PetscTruth*)&i);
1029   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IS is not sorted");
1030 
1031   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
1032   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
1033   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
1034 
1035   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
1036     /* special case of contiguous rows */
1037     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
1038     starts = lens + ncols;
1039     /* loop over new rows determining lens and starting points */
1040     for (i=0; i<nrows; i++) {
1041       kstart  = ai[irow[i]]+shift;
1042       kend    = kstart + ailen[irow[i]];
1043       for ( k=kstart; k<kend; k++ ) {
1044         if (aj[k]+shift >= first) {
1045           starts[i] = k;
1046           break;
1047 	}
1048       }
1049       sum = 0;
1050       while (k < kend) {
1051         if (aj[k++]+shift >= first+ncols) break;
1052         sum++;
1053       }
1054       lens[i] = sum;
1055     }
1056     /* create submatrix */
1057     if (scall == MAT_REUSE_MATRIX) {
1058       int n_cols,n_rows;
1059       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1060       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1061       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
1062       C = *B;
1063     }
1064     else {
1065       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
1066     }
1067     c = (Mat_SeqAIJ*) C->data;
1068 
1069     /* loop over rows inserting into submatrix */
1070     a_new    = c->a;
1071     j_new    = c->j;
1072     i_new    = c->i;
1073     i_new[0] = -shift;
1074     for (i=0; i<nrows; i++) {
1075       ii    = starts[i];
1076       lensi = lens[i];
1077       for ( k=0; k<lensi; k++ ) {
1078         *j_new++ = aj[ii+k] - first;
1079       }
1080       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1081       a_new      += lensi;
1082       i_new[i+1]  = i_new[i] + lensi;
1083       c->ilen[i]  = lensi;
1084     }
1085     PetscFree(lens);
1086   }
1087   else {
1088     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
1089     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
1090     ssmap = smap + shift;
1091     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1092     PetscMemzero(smap,oldcols*sizeof(int));
1093     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
1094     /* determine lens of each row */
1095     for (i=0; i<nrows; i++) {
1096       kstart  = ai[irow[i]]+shift;
1097       kend    = kstart + a->ilen[irow[i]];
1098       lens[i] = 0;
1099       for ( k=kstart; k<kend; k++ ) {
1100         if (ssmap[aj[k]]) {
1101           lens[i]++;
1102         }
1103       }
1104     }
1105     /* Create and fill new matrix */
1106     if (scall == MAT_REUSE_MATRIX) {
1107       c = (Mat_SeqAIJ *)((*B)->data);
1108 
1109       if (c->m  != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1110       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1111         SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros");
1112       }
1113       PetscMemzero(c->ilen,c->m*sizeof(int));
1114       C = *B;
1115     }
1116     else {
1117       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
1118     }
1119     c = (Mat_SeqAIJ *)(C->data);
1120     for (i=0; i<nrows; i++) {
1121       row    = irow[i];
1122       nznew  = 0;
1123       kstart = ai[row]+shift;
1124       kend   = kstart + a->ilen[row];
1125       mat_i  = c->i[i]+shift;
1126       mat_j  = c->j + mat_i;
1127       mat_a  = c->a + mat_i;
1128       mat_ilen = c->ilen + i;
1129       for ( k=kstart; k<kend; k++ ) {
1130         if ((tcol=ssmap[a->j[k]])) {
1131           *mat_j++ = tcol - (!shift);
1132           *mat_a++ = a->a[k];
1133           (*mat_ilen)++;
1134 
1135         }
1136       }
1137     }
1138     /* Free work space */
1139     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1140     PetscFree(smap); PetscFree(lens);
1141   }
1142   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1143   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1144 
1145   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1146   *B = C;
1147   return 0;
1148 }
1149 
1150 /*
1151      note: This can only work for identity for row and col. It would
1152    be good to check this and otherwise generate an error.
1153 */
1154 static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1155 {
1156   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1157   int        ierr;
1158   Mat        outA;
1159 
1160   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1161 
1162   outA          = inA;
1163   inA->factor   = FACTOR_LU;
1164   a->row        = row;
1165   a->col        = col;
1166 
1167   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
1168 
1169   if (!a->diag) {
1170     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
1171   }
1172   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1173   return 0;
1174 }
1175 
1176 #include "pinclude/plapack.h"
1177 static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1178 {
1179   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1180   int        one = 1;
1181   BLscal_( &a->nz, alpha, a->a, &one );
1182   PLogFlops(a->nz);
1183   return 0;
1184 }
1185 
1186 static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1187                                     Mat **B)
1188 {
1189   int ierr,i;
1190 
1191   if (scall == MAT_INITIAL_MATRIX) {
1192     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1193   }
1194 
1195   for ( i=0; i<n; i++ ) {
1196     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1197   }
1198   return 0;
1199 }
1200 
1201 static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
1202 {
1203   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1204   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
1205   int        start, end, *ai, *aj;
1206   char       *table;
1207   shift = a->indexshift;
1208   m     = a->m;
1209   ai    = a->i;
1210   aj    = a->j+shift;
1211 
1212   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
1213 
1214   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
1215   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1216 
1217   for ( i=0; i<is_max; i++ ) {
1218     /* Initialise the two local arrays */
1219     isz  = 0;
1220     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1221 
1222                 /* Extract the indices, assume there can be duplicate entries */
1223     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
1224     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1225 
1226     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1227     for ( j=0; j<n ; ++j){
1228       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
1229     }
1230     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
1231     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1232 
1233     k = 0;
1234     for ( j=0; j<ov; j++){ /* for each overlap*/
1235       n = isz;
1236       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1237         row   = nidx[k];
1238         start = ai[row];
1239         end   = ai[row+1];
1240         for ( l = start; l<end ; l++){
1241           val = aj[l] + shift;
1242           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1243         }
1244       }
1245     }
1246     ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1247   }
1248   PetscFree(table);
1249   PetscFree(nidx);
1250   return 0;
1251 }
1252 
1253 int MatPrintHelp_SeqAIJ(Mat A)
1254 {
1255   static int called = 0;
1256   MPI_Comm   comm = A->comm;
1257 
1258   if (called) return 0; else called = 1;
1259   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1260   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
1261   PetscPrintf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
1262   PetscPrintf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
1263   PetscPrintf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1264 #if defined(HAVE_ESSL)
1265   PetscPrintf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1266 #endif
1267   return 0;
1268 }
1269 static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1270 /* -------------------------------------------------------------------*/
1271 static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
1272        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1273        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1274        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
1275        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
1276        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
1277        MatLUFactor_SeqAIJ,0,
1278        MatRelax_SeqAIJ,
1279        MatTranspose_SeqAIJ,
1280        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1281        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
1282        0,MatAssemblyEnd_SeqAIJ,
1283        MatCompress_SeqAIJ,
1284        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
1285        MatGetReordering_SeqAIJ,
1286        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
1287        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
1288        MatILUFactorSymbolic_SeqAIJ,0,
1289        0,0,MatConvert_SeqAIJ,
1290        MatGetSubMatrix_SeqAIJ,0,
1291        MatConvertSameType_SeqAIJ,0,0,
1292        MatILUFactor_SeqAIJ,0,0,
1293        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1294        MatGetValues_SeqAIJ,0,
1295        MatPrintHelp_SeqAIJ,
1296        MatScale_SeqAIJ};
1297 
1298 extern int MatUseSuperLU_SeqAIJ(Mat);
1299 extern int MatUseEssl_SeqAIJ(Mat);
1300 extern int MatUseDXML_SeqAIJ(Mat);
1301 
1302 /*@C
1303    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
1304    (the default parallel PETSc format).  For good matrix assembly performance
1305    the user should preallocate the matrix storage by setting the parameter nz
1306    (or nzz).  By setting these parameters accurately, performance can be
1307    increased by more than a factor of 50.
1308 
1309    Input Parameters:
1310 .  comm - MPI communicator, set to MPI_COMM_SELF
1311 .  m - number of rows
1312 .  n - number of columns
1313 .  nz - number of nonzeros per row (same for all rows)
1314 .  nzz - number of nonzeros per row or null (possibly different for each row)
1315 
1316    Output Parameter:
1317 .  A - the matrix
1318 
1319    Notes:
1320    The AIJ format (also called the Yale sparse matrix format or
1321    compressed row storage), is fully compatible with standard Fortran 77
1322    storage.  That is, the stored row and column indices can begin at
1323    either one (as in Fortran) or zero.  See the users manual for details.
1324 
1325    Specify the preallocated storage with either nz or nnz (not both).
1326    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1327    allocation.  For additional details, see the users manual chapter on
1328    matrices and the file $(PETSC_DIR)/Performance.
1329 
1330    By default, this format uses inodes (identical nodes) when possible, to
1331    improve numerical efficiency of Matrix vector products and solves. We
1332    search for consecutive rows with the same nonzero structure, thereby
1333    reusing matrix information to achieve increased efficiency.
1334 
1335    Options Database Keys:
1336 $    -mat_aij_no_inode  - Do not use inodes
1337 $    -mat_aij_inode_limit <limit> - Set inode limit.
1338 $        (max limit=5)
1339 $    -mat_aij_oneindex - Internally use indexing starting at 1
1340 $        rather than 0.  Note: When calling MatSetValues(),
1341 $        the user still MUST index entries starting at 0!
1342 
1343 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
1344 @*/
1345 int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
1346 {
1347   Mat        B;
1348   Mat_SeqAIJ *b;
1349   int        i,len,ierr, flg;
1350 
1351   *A      = 0;
1352   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1353   PLogObjectCreate(B);
1354   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1355   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1356   B->destroy          = MatDestroy_SeqAIJ;
1357   B->view             = MatView_SeqAIJ;
1358   B->factor           = 0;
1359   B->lupivotthreshold = 1.0;
1360   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
1361                           &flg); CHKERRQ(ierr);
1362   b->ilu_preserve_row_sums = PETSC_FALSE;
1363   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
1364                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1365   b->row              = 0;
1366   b->col              = 0;
1367   b->indexshift       = 0;
1368   b->reallocs         = 0;
1369   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
1370   if (flg) b->indexshift = -1;
1371 
1372   b->m       = m;
1373   b->n       = n;
1374   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1375   if (nnz == PETSC_NULL) {
1376     if (nz == PETSC_DEFAULT) nz = 10;
1377     else if (nz <= 0)        nz = 1;
1378     for ( i=0; i<m; i++ ) b->imax[i] = nz;
1379     nz = nz*m;
1380   }
1381   else {
1382     nz = 0;
1383     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
1384   }
1385 
1386   /* allocate the matrix space */
1387   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
1388   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1389   b->j  = (int *) (b->a + nz);
1390   PetscMemzero(b->j,nz*sizeof(int));
1391   b->i  = b->j + nz;
1392   b->singlemalloc = 1;
1393 
1394   b->i[0] = -b->indexshift;
1395   for (i=1; i<m+1; i++) {
1396     b->i[i] = b->i[i-1] + b->imax[i-1];
1397   }
1398 
1399   /* b->ilen will count nonzeros in each row so far. */
1400   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1401   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1402   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
1403 
1404   b->nz               = 0;
1405   b->maxnz            = nz;
1406   b->sorted           = 0;
1407   b->roworiented      = 1;
1408   b->nonew            = 0;
1409   b->diag             = 0;
1410   b->solve_work       = 0;
1411   b->spptr            = 0;
1412   b->inode.node_count = 0;
1413   b->inode.size       = 0;
1414   b->inode.limit      = 5;
1415   b->inode.max_limit  = 5;
1416 
1417   *A = B;
1418   /*  SuperLU is not currently supported through PETSc */
1419 #if defined(HAVE_SUPERLU)
1420   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
1421   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
1422 #endif
1423   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
1424   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
1425   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
1426   if (flg) {
1427     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1428     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
1429   }
1430   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1431   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1432   return 0;
1433 }
1434 
1435 int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
1436 {
1437   Mat        C;
1438   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1439   int        i,len, m = a->m,shift = a->indexshift;
1440 
1441   *B = 0;
1442   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1443   PLogObjectCreate(C);
1444   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1445   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1446   C->destroy    = MatDestroy_SeqAIJ;
1447   C->view       = MatView_SeqAIJ;
1448   C->factor     = A->factor;
1449   c->row        = 0;
1450   c->col        = 0;
1451   c->indexshift = shift;
1452   C->assembled  = PETSC_TRUE;
1453 
1454   c->m          = a->m;
1455   c->n          = a->n;
1456 
1457   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1458   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
1459   for ( i=0; i<m; i++ ) {
1460     c->imax[i] = a->imax[i];
1461     c->ilen[i] = a->ilen[i];
1462   }
1463 
1464   /* allocate the matrix space */
1465   c->singlemalloc = 1;
1466   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1467   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1468   c->j  = (int *) (c->a + a->i[m] + shift);
1469   c->i  = c->j + a->i[m] + shift;
1470   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
1471   if (m > 0) {
1472     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1473     if (cpvalues == COPY_VALUES) {
1474       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
1475     }
1476   }
1477 
1478   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1479   c->sorted      = a->sorted;
1480   c->roworiented = a->roworiented;
1481   c->nonew       = a->nonew;
1482   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
1483 
1484   if (a->diag) {
1485     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1486     PLogObjectMemory(C,(m+1)*sizeof(int));
1487     for ( i=0; i<m; i++ ) {
1488       c->diag[i] = a->diag[i];
1489     }
1490   }
1491   else c->diag          = 0;
1492   c->inode.limit        = a->inode.limit;
1493   c->inode.max_limit    = a->inode.max_limit;
1494   if (a->inode.size){
1495     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1496     c->inode.node_count = a->inode.node_count;
1497     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1498   } else {
1499     c->inode.size       = 0;
1500     c->inode.node_count = 0;
1501   }
1502   c->nz                 = a->nz;
1503   c->maxnz              = a->maxnz;
1504   c->solve_work         = 0;
1505   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1506 
1507   *B = C;
1508   return 0;
1509 }
1510 
1511 int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
1512 {
1513   Mat_SeqAIJ   *a;
1514   Mat          B;
1515   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1516   MPI_Comm     comm;
1517 
1518   PetscObjectGetComm((PetscObject) viewer,&comm);
1519   MPI_Comm_size(comm,&size);
1520   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
1521   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
1522   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1523   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
1524   M = header[1]; N = header[2]; nz = header[3];
1525 
1526   /* read in row lengths */
1527   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1528   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
1529 
1530   /* create our matrix */
1531   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1532   B = *A;
1533   a = (Mat_SeqAIJ *) B->data;
1534   shift = a->indexshift;
1535 
1536   /* read in column indices and adjust for Fortran indexing*/
1537   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
1538   if (shift) {
1539     for ( i=0; i<nz; i++ ) {
1540       a->j[i] += 1;
1541     }
1542   }
1543 
1544   /* read in nonzero values */
1545   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
1546 
1547   /* set matrix "i" values */
1548   a->i[0] = -shift;
1549   for ( i=1; i<= M; i++ ) {
1550     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1551     a->ilen[i-1] = rowlengths[i-1];
1552   }
1553   PetscFree(rowlengths);
1554 
1555   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1556   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1557   return 0;
1558 }
1559 
1560 static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
1561 {
1562   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
1563 
1564   if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type");
1565 
1566   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
1567   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1568       (a->indexshift != b->indexshift)) {
1569     *flg = PETSC_FALSE; return 0;
1570   }
1571 
1572   /* if the a->i are the same */
1573   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
1574     *flg = PETSC_FALSE; return 0;
1575   }
1576 
1577   /* if a->j are the same */
1578   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
1579     *flg = PETSC_FALSE; return 0;
1580   }
1581 
1582   /* if a->a are the same */
1583   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
1584     *flg = PETSC_FALSE; return 0;
1585   }
1586   *flg = PETSC_TRUE;
1587   return 0;
1588 
1589 }
1590