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