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