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