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