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