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