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