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