xref: /petsc/src/mat/impls/aij/seq/aij.c (revision ddd1276785a76ea4a3b2575e2a65f0ba44dae580)
1 #ifndef lint
2 static char vcid[] = "$Id: aij.c,v 1.97 1995/10/11 22:09:34 curfman 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 == COLUMNS_SORTED)            a->sorted      = 1;
330   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
331   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
332   else if (op == ROWS_SORTED ||
333            op == SYMMETRIC_MATRIX ||
334            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
335            op == YES_NEW_DIAGONALS)
336     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
337   else if (op == COLUMN_ORIENTED)
338     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED not supported");}
339   else if (op == NO_NEW_DIAGONALS)
340     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS not supported");}
341   else
342     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:Option not supported");}
343   return 0;
344 }
345 
346 static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
347 {
348   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
349   int        i,j, n,shift = a->indexshift;
350   Scalar     *x, zero = 0.0;
351 
352   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
353   VecSet(&zero,v);
354   VecGetArray(v,&x); VecGetLocalSize(v,&n);
355   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
356   for ( i=0; i<a->m; i++ ) {
357     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
358       if (a->j[j]+shift == i) {
359         x[i] = a->a[j];
360         break;
361       }
362     }
363   }
364   return 0;
365 }
366 
367 /* -------------------------------------------------------*/
368 /* Should check that shapes of vectors and matrices match */
369 /* -------------------------------------------------------*/
370 static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
371 {
372   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
373   Scalar     *x, *y, *v, alpha;
374   int        m = a->m, n, i, *idx, shift = a->indexshift;
375 
376   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
377   VecGetArray(xx,&x); VecGetArray(yy,&y);
378   PetscZero(y,a->n*sizeof(Scalar));
379   y = y + shift; /* shift for Fortran start by 1 indexing */
380   for ( i=0; i<m; i++ ) {
381     idx   = a->j + a->i[i] + shift;
382     v     = a->a + a->i[i] + shift;
383     n     = a->i[i+1] - a->i[i];
384     alpha = x[i];
385     while (n-->0) {y[*idx++] += alpha * *v++;}
386   }
387   PLogFlops(2*a->nz - a->n);
388   return 0;
389 }
390 static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
391 {
392   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
393   Scalar     *x, *y, *v, alpha;
394   int        m = a->m, n, i, *idx,shift = a->indexshift;
395 
396   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
397   VecGetArray(xx,&x); VecGetArray(yy,&y);
398   if (zz != yy) VecCopy(zz,yy);
399   y = y + shift; /* shift for Fortran start by 1 indexing */
400   for ( i=0; i<m; i++ ) {
401     idx   = a->j + a->i[i] + shift;
402     v     = a->a + a->i[i] + shift;
403     n     = a->i[i+1] - a->i[i];
404     alpha = x[i];
405     while (n-->0) {y[*idx++] += alpha * *v++;}
406   }
407   return 0;
408 }
409 
410 static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
411 {
412   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
413   Scalar     *x, *y, *v, sum;
414   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
415 
416   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
417   VecGetArray(xx,&x); VecGetArray(yy,&y);
418   x = x + shift; /* shift for Fortran start by 1 indexing */
419   idx  = a->j;
420   v    = a->a;
421   ii   = a->i;
422 #if defined(PARCH_rs6000)
423 #pragma disjoint (*x,*v,*y)
424 #endif
425   for ( i=0; i<m; i++ ) {
426     n    = ii[1] - ii[0]; ii++;
427     sum  = 0.0;
428     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
429     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
430     while (n--) sum += *v++ * x[*idx++];
431     y[i] = sum;
432   }
433   PLogFlops(2*a->nz - m);
434   return 0;
435 }
436 
437 static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
438 {
439   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
440   Scalar     *x, *y, *z, *v, sum;
441   int        m = a->m, n, i, *idx, shift = a->indexshift;
442 
443   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
444   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
445   x = x + shift; /* shift for Fortran start by 1 indexing */
446   for ( i=0; i<m; i++ ) {
447     idx  = a->j + a->i[i] + shift;
448     v    = a->a + a->i[i] + shift;
449     n    = a->i[i+1] - a->i[i];
450     sum  = y[i];
451     SPARSEDENSEDOT(sum,x,v,idx,n);
452     z[i] = sum;
453   }
454   PLogFlops(2*a->nz);
455   return 0;
456 }
457 
458 /*
459      Adds diagonal pointers to sparse matrix structure.
460 */
461 
462 int MatMarkDiag_SeqAIJ(Mat A)
463 {
464   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
465   int        i,j, *diag, m = a->m,shift = a->indexshift;
466 
467   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
468   diag = (int *) PETSCMALLOC( (m+1)*sizeof(int)); CHKPTRQ(diag);
469   PLogObjectMemory(A,(m+1)*sizeof(int));
470   for ( i=0; i<a->m; i++ ) {
471     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
472       if (a->j[j]+shift == i) {
473         diag[i] = j - shift;
474         break;
475       }
476     }
477   }
478   a->diag = diag;
479   return 0;
480 }
481 
482 static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
483                            double fshift,int its,Vec xx)
484 {
485   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
486   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
487   int        ierr, *idx, *diag,n = a->n, m = a->m, i;
488   int        shift = a->indexshift;
489 
490   VecGetArray(xx,&x); VecGetArray(bb,&b);
491   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
492   diag = a->diag;
493   xs   = x + shift; /* shifted by one for index start of a or a->j*/
494   if (flag == SOR_APPLY_UPPER) {
495    /* apply ( U + D/omega) to the vector */
496     bs = b + shift;
497     for ( i=0; i<m; i++ ) {
498         d    = fshift + a->a[diag[i] + shift];
499         n    = a->i[i+1] - diag[i] - 1;
500         idx  = a->j + diag[i] + (!shift);
501         v    = a->a + diag[i] + (!shift);
502         sum  = b[i]*d/omega;
503         SPARSEDENSEDOT(sum,bs,v,idx,n);
504         x[i] = sum;
505     }
506     return 0;
507   }
508   if (flag == SOR_APPLY_LOWER) {
509     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
510   }
511   else if (flag & SOR_EISENSTAT) {
512     /* Let  A = L + U + D; where L is lower trianglar,
513     U is upper triangular, E is diagonal; This routine applies
514 
515             (L + E)^{-1} A (U + E)^{-1}
516 
517     to a vector efficiently using Eisenstat's trick. This is for
518     the case of SSOR preconditioner, so E is D/omega where omega
519     is the relaxation factor.
520     */
521     t = (Scalar *) PETSCMALLOC( m*sizeof(Scalar) ); CHKPTRQ(t);
522     scale = (2.0/omega) - 1.0;
523 
524     /*  x = (E + U)^{-1} b */
525     for ( i=m-1; i>=0; i-- ) {
526       d    = fshift + a->a[diag[i] + shift];
527       n    = a->i[i+1] - diag[i] - 1;
528       idx  = a->j + diag[i] + (!shift);
529       v    = a->a + diag[i] + (!shift);
530       sum  = b[i];
531       SPARSEDENSEMDOT(sum,xs,v,idx,n);
532       x[i] = omega*(sum/d);
533     }
534 
535     /*  t = b - (2*E - D)x */
536     v = a->a;
537     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
538 
539     /*  t = (E + L)^{-1}t */
540     ts = t + shift; /* shifted by one for index start of a or a->j*/
541     diag = a->diag;
542     for ( i=0; i<m; i++ ) {
543       d    = fshift + a->a[diag[i]+shift];
544       n    = diag[i] - a->i[i];
545       idx  = a->j + a->i[i] + shift;
546       v    = a->a + a->i[i] + shift;
547       sum  = t[i];
548       SPARSEDENSEMDOT(sum,ts,v,idx,n);
549       t[i] = omega*(sum/d);
550     }
551 
552     /*  x = x + t */
553     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
554     PETSCFREE(t);
555     return 0;
556   }
557   if (flag & SOR_ZERO_INITIAL_GUESS) {
558     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
559       for ( i=0; i<m; i++ ) {
560         d    = fshift + a->a[diag[i]+shift];
561         n    = diag[i] - a->i[i];
562         idx  = a->j + a->i[i] + shift;
563         v    = a->a + a->i[i] + shift;
564         sum  = b[i];
565         SPARSEDENSEMDOT(sum,xs,v,idx,n);
566         x[i] = omega*(sum/d);
567       }
568       xb = x;
569     }
570     else xb = b;
571     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
572         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
573       for ( i=0; i<m; i++ ) {
574         x[i] *= a->a[diag[i]+shift];
575       }
576     }
577     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
578       for ( i=m-1; i>=0; i-- ) {
579         d    = fshift + a->a[diag[i] + shift];
580         n    = a->i[i+1] - diag[i] - 1;
581         idx  = a->j + diag[i] + (!shift);
582         v    = a->a + diag[i] + (!shift);
583         sum  = xb[i];
584         SPARSEDENSEMDOT(sum,xs,v,idx,n);
585         x[i] = omega*(sum/d);
586       }
587     }
588     its--;
589   }
590   while (its--) {
591     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
592       for ( i=0; i<m; i++ ) {
593         d    = fshift + a->a[diag[i]+shift];
594         n    = a->i[i+1] - a->i[i];
595         idx  = a->j + a->i[i] + shift;
596         v    = a->a + a->i[i] + shift;
597         sum  = b[i];
598         SPARSEDENSEMDOT(sum,xs,v,idx,n);
599         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
600       }
601     }
602     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
603       for ( i=m-1; i>=0; i-- ) {
604         d    = fshift + a->a[diag[i] + shift];
605         n    = a->i[i+1] - a->i[i];
606         idx  = a->j + a->i[i] + shift;
607         v    = a->a + a->i[i] + shift;
608         sum  = b[i];
609         SPARSEDENSEMDOT(sum,xs,v,idx,n);
610         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
611       }
612     }
613   }
614   return 0;
615 }
616 
617 static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,
618                              int *nzalloc,int *mem)
619 {
620   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
621   *nz      = a->nz;
622   *nzalloc = a->maxnz;
623   *mem     = (int)A->mem;
624   return 0;
625 }
626 
627 extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
628 extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
629 extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
630 extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
631 extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
632 extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
633 extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
634 
635 static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
636 {
637   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
638   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
639 
640   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
641   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
642   if (diag) {
643     for ( i=0; i<N; i++ ) {
644       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
645       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
646         a->ilen[rows[i]] = 1;
647         a->a[a->i[rows[i]]+shift] = *diag;
648         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
649       }
650       else {
651         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
652         CHKERRQ(ierr);
653       }
654     }
655   }
656   else {
657     for ( i=0; i<N; i++ ) {
658       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
659       a->ilen[rows[i]] = 0;
660     }
661   }
662   ISRestoreIndices(is,&rows);
663   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
664   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
665   return 0;
666 }
667 
668 static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
669 {
670   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
671   *m = a->m; *n = a->n;
672   return 0;
673 }
674 
675 static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
676 {
677   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
678   *m = 0; *n = a->m;
679   return 0;
680 }
681 static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
682 {
683   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
684   int        *itmp,i,ierr,shift = a->indexshift;
685 
686   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
687 
688   if (!a->assembled) {
689     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
690     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
691   }
692   *nz = a->i[row+1] - a->i[row];
693   if (v) *v = a->a + a->i[row] + shift;
694   if (idx) {
695     if (*nz) {
696       itmp = a->j + a->i[row] + shift;
697       *idx = (int *) PETSCMALLOC( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
698       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
699     }
700     else *idx = 0;
701   }
702   return 0;
703 }
704 
705 static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
706 {
707   if (idx) {if (*idx) PETSCFREE(*idx);}
708   return 0;
709 }
710 
711 static int MatNorm_SeqAIJ(Mat A,MatNormType type,double *norm)
712 {
713   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
714   Scalar     *v = a->a;
715   double     sum = 0.0;
716   int        i, j,shift = a->indexshift;
717 
718   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
719   if (type == NORM_FROBENIUS) {
720     for (i=0; i<a->nz; i++ ) {
721 #if defined(PETSC_COMPLEX)
722       sum += real(conj(*v)*(*v)); v++;
723 #else
724       sum += (*v)*(*v); v++;
725 #endif
726     }
727     *norm = sqrt(sum);
728   }
729   else if (type == NORM_1) {
730     double *tmp;
731     int    *jj = a->j;
732     tmp = (double *) PETSCMALLOC( a->n*sizeof(double) ); CHKPTRQ(tmp);
733     PetscZero(tmp,a->n*sizeof(double));
734     *norm = 0.0;
735     for ( j=0; j<a->nz; j++ ) {
736 #if defined(PETSC_COMPLEX)
737         tmp[*jj++ + shift] += abs(*v++);
738 #else
739         tmp[*jj++ + shift] += fabs(*v++);
740 #endif
741     }
742     for ( j=0; j<a->n; j++ ) {
743       if (tmp[j] > *norm) *norm = tmp[j];
744     }
745     PETSCFREE(tmp);
746   }
747   else if (type == NORM_INFINITY) {
748     *norm = 0.0;
749     for ( j=0; j<a->m; j++ ) {
750       v = a->a + a->i[j] + shift;
751       sum = 0.0;
752       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
753 #if defined(PETSC_COMPLEX)
754         sum += abs(*v); v++;
755 #else
756         sum += fabs(*v); v++;
757 #endif
758       }
759       if (sum > *norm) *norm = sum;
760     }
761   }
762   else {
763     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
764   }
765   return 0;
766 }
767 
768 static int MatTranspose_SeqAIJ(Mat A,Mat *B)
769 {
770   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
771   Mat        C;
772   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
773   Scalar     *array = a->a;
774   int        shift = a->indexshift;
775 
776   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
777   col = (int *) PETSCMALLOC((1+a->n)*sizeof(int)); CHKPTRQ(col);
778   PetscZero(col,(1+a->n)*sizeof(int));
779   if (shift) {
780     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
781   }
782   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
783   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
784   PETSCFREE(col);
785   for ( i=0; i<m; i++ ) {
786     len = ai[i+1]-ai[i];
787     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
788     array += len; aj += len;
789   }
790   if (shift) {
791     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
792   }
793 
794   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
795   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
796 
797   if (B) {
798     *B = C;
799   } else {
800     /* This isn't really an in-place transpose */
801     PETSCFREE(a->a);
802     if (!a->singlemalloc) {PETSCFREE(a->i); PETSCFREE(a->j);}
803     if (a->diag) PETSCFREE(a->diag);
804     if (a->ilen) PETSCFREE(a->ilen);
805     if (a->imax) PETSCFREE(a->imax);
806     if (a->solve_work) PETSCFREE(a->solve_work);
807     PETSCFREE(a);
808     PetscMemcpy(A,C,sizeof(struct _Mat));
809     PETSCHEADERDESTROY(C);
810   }
811   return 0;
812 }
813 
814 static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
815 {
816   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
817   Scalar     *l,*r,x,*v;
818   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj;
819   int        shift = a->indexshift;
820 
821   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
822   if (ll) {
823     VecGetArray(ll,&l); VecGetSize(ll,&m);
824     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
825     v = a->a;
826     for ( i=0; i<m; i++ ) {
827       x = l[i];
828       M = a->i[i+1] - a->i[i];
829       for ( j=0; j<M; j++ ) { (*v++) *= x;}
830     }
831   }
832   if (rr) {
833     VecGetArray(rr,&r); VecGetSize(rr,&n);
834     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
835     v = a->a; jj = a->j;
836     for ( i=0; i<nz; i++ ) {
837       (*v++) *= r[*jj++ + shift];
838     }
839   }
840   return 0;
841 }
842 
843 static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,Mat *B)
844 {
845   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c;
846   int        nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
847   int        *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
848   int        first,step,*starts,*j_new,*i_new;
849   Scalar     *vwork,*a_new;
850   Mat        C;
851 
852   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
853   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
854   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
855   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
856 
857   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
858     /* special case of contiguous rows */
859     lens   = (int *) PETSCMALLOC((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
860     starts = lens + ncols;
861     /* loop over new rows determining lens and starting points */
862     for (i=0; i<nrows; i++) {
863       kstart  = a->i[irow[i]]+shift;
864       kend    = kstart + a->ilen[irow[i]];
865       for ( k=kstart; k<kend; k++ ) {
866         if (a->j[k] >= first) {
867           starts[i] = k;
868           break;
869 	}
870       }
871       lens[i] = 0;
872       while (k < kend) {
873         if (a->j[k++] >= first+ncols) break;
874         lens[i]++;
875       }
876     }
877     /* create submatrix */
878     ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
879     c = (Mat_SeqAIJ*) C->data;
880 
881     /* loop over rows inserting into submatrix */
882     a_new    = c->a;
883     j_new    = c->j;
884     i_new    = c->i;
885     i_new[0] = -shift;
886     for (i=0; i<nrows; i++) {
887       for ( k=0; k<lens[i]; k++ ) {
888         *j_new++ = a->j[k+starts[i]] - first;
889       }
890       PetscMemcpy(a_new,a->a + starts[i],lens[i]*sizeof(Scalar));
891       a_new      += lens[i];
892       i_new[i+1]  = i_new[i] + lens[i];
893       c->ilen[i]  = lens[i];
894     }
895     PETSCFREE(lens);
896   }
897   else {
898     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
899     smap  = (int *) PETSCMALLOC((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
900     ssmap = smap + shift;
901     cwork = (int *) PETSCMALLOC((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork);
902     lens  = cwork + ncols;
903     vwork = (Scalar *) PETSCMALLOC((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
904     PetscZero(smap,oldcols*sizeof(int));
905     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
906     /* determine lens of each row */
907     for (i=0; i<nrows; i++) {
908       kstart  = a->i[irow[i]]+shift;
909       kend    = kstart + a->ilen[irow[i]];
910       lens[i] = 0;
911       for ( k=kstart; k<kend; k++ ) {
912         if (ssmap[a->j[k]]) {
913           lens[i]++;
914         }
915       }
916     }
917     /* Create and fill new matrix */
918     ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
919     for (i=0; i<nrows; i++) {
920       nznew  = 0;
921       kstart = a->i[irow[i]]+shift;
922       kend   = kstart + a->ilen[irow[i]];
923       for ( k=kstart; k<kend; k++ ) {
924         if (ssmap[a->j[k]]) {
925           cwork[nznew]   = ssmap[a->j[k]] - 1;
926           vwork[nznew++] = a->a[k];
927         }
928       }
929       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
930     }
931     /* Free work space */
932     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
933     PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
934   }
935   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
936   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
937 
938   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
939   *B = C;
940   return 0;
941 }
942 
943 /* -------------------------------------------------------------------*/
944 extern int MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,int,Mat *);
945 extern int MatConvert_SeqAIJ(Mat,MatType,Mat *);
946 static int MatCopyPrivate_SeqAIJ(Mat,Mat *);
947 
948 static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
949        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
950        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
951        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
952        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
953        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
954        MatLUFactor_SeqAIJ,0,
955        MatRelax_SeqAIJ,
956        MatTranspose_SeqAIJ,
957        MatGetInfo_SeqAIJ,0,
958        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
959        0,MatAssemblyEnd_SeqAIJ,
960        MatCompress_SeqAIJ,
961        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
962        MatGetReordering_SeqAIJ,
963        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
964        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
965        MatILUFactorSymbolic_SeqAIJ,0,
966        0,0,MatConvert_SeqAIJ,
967        MatGetSubMatrix_SeqAIJ,0,
968        MatCopyPrivate_SeqAIJ};
969 
970 extern int MatUseSuperLU_SeqAIJ(Mat);
971 extern int MatUseEssl_SeqAIJ(Mat);
972 extern int MatUseDXML_SeqAIJ(Mat);
973 
974 /*@C
975    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
976    (the default uniprocessor PETSc format).
977 
978    Input Parameters:
979 .  comm - MPI communicator, set to MPI_COMM_SELF
980 .  m - number of rows
981 .  n - number of columns
982 .  nz - number of nonzeros per row (same for all rows)
983 .  nzz - number of nonzeros per row or null (possibly different for each row)
984 
985    Output Parameter:
986 .  A - the matrix
987 
988    Notes:
989    The AIJ format (also called the Yale sparse matrix format or
990    compressed row storage), is fully compatible with standard Fortran 77
991    storage.  That is, the stored row and column indices begin at
992    one, not zero.
993 
994    Specify the preallocated storage with either nz or nnz (not both).
995    Set both nz and nnz to zero for PETSc to control dynamic memory
996    allocation.
997 
998 .keywords: matrix, aij, compressed row, sparse
999 
1000 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
1001 @*/
1002 int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
1003 {
1004   Mat        B;
1005   Mat_SeqAIJ *b;
1006   int        i,len,ierr;
1007   *A      = 0;
1008   PETSCHEADERCREATE(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1009   PLogObjectCreate(B);
1010   B->data               = (void *) (b = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(b);
1011   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1012   B->destroy          = MatDestroy_SeqAIJ;
1013   B->view             = MatView_SeqAIJ;
1014   B->factor           = 0;
1015   B->lupivotthreshold = 1.0;
1016   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1017   b->row              = 0;
1018   b->col              = 0;
1019   b->indexshift       = 0;
1020   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
1021 
1022   b->m       = m;
1023   b->n       = n;
1024   b->imax    = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1025   if (!nnz) {
1026     if (nz <= 0) nz = 1;
1027     for ( i=0; i<m; i++ ) b->imax[i] = nz;
1028     nz = nz*m;
1029   }
1030   else {
1031     nz = 0;
1032     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
1033   }
1034 
1035   /* allocate the matrix space */
1036   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
1037   b->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(b->a);
1038   b->j  = (int *) (b->a + nz);
1039   PetscZero(b->j,nz*sizeof(int));
1040   b->i  = b->j + nz;
1041   b->singlemalloc = 1;
1042 
1043   b->i[0] = -b->indexshift;
1044   for (i=1; i<m+1; i++) {
1045     b->i[i] = b->i[i-1] + b->imax[i-1];
1046   }
1047 
1048   /* b->ilen will count nonzeros in each row so far. */
1049   b->ilen = (int *) PETSCMALLOC((m+1)*sizeof(int));
1050   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1051   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
1052 
1053   b->nz          = 0;
1054   b->maxnz       = nz;
1055   b->sorted      = 0;
1056   b->roworiented = 1;
1057   b->nonew       = 0;
1058   b->diag        = 0;
1059   b->assembled   = 0;
1060   b->solve_work  = 0;
1061 
1062   *A = B;
1063   if (OptionsHasName(0,"-mat_aij_superlu")) {
1064     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
1065   }
1066   if (OptionsHasName(0,"-mat_aij_essl")) {
1067     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
1068   }
1069   if (OptionsHasName(0,"-mat_aij_dxml")) {
1070     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1071     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
1072   }
1073 
1074   return 0;
1075 }
1076 
1077 static int MatCopyPrivate_SeqAIJ(Mat A,Mat *B)
1078 {
1079   Mat        C;
1080   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1081   int        i,len, m = a->m;
1082   int        shift = a->indexshift;
1083   *B      = 0;
1084 
1085   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
1086   PETSCHEADERCREATE(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1087   PLogObjectCreate(C);
1088   C->data       = (void *) (c = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(c);
1089   PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps));
1090   C->destroy    = MatDestroy_SeqAIJ;
1091   C->view       = MatView_SeqAIJ;
1092   C->factor     = A->factor;
1093   c->row        = 0;
1094   c->col        = 0;
1095   c->indexshift = shift;
1096 
1097   c->m          = a->m;
1098   c->n          = a->n;
1099 
1100   c->imax       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1101   c->ilen       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
1102   for ( i=0; i<m; i++ ) {
1103     c->imax[i] = a->imax[i];
1104     c->ilen[i] = a->ilen[i];
1105   }
1106 
1107   /* allocate the matrix space */
1108   c->singlemalloc = 1;
1109   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1110   c->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(c->a);
1111   c->j  = (int *) (c->a + a->i[m] + shift);
1112   c->i  = c->j + a->i[m] + shift;
1113   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
1114   if (m > 0) {
1115     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1116     PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
1117   }
1118 
1119   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1120   c->sorted      = a->sorted;
1121   c->roworiented = a->roworiented;
1122   c->nonew       = a->nonew;
1123 
1124   if (a->diag) {
1125     c->diag = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1126     PLogObjectMemory(C,(m+1)*sizeof(int));
1127     for ( i=0; i<m; i++ ) {
1128       c->diag[i] = a->diag[i];
1129     }
1130   }
1131   else c->diag        = 0;
1132   c->assembled        = 1;
1133   c->nz               = a->nz;
1134   c->maxnz            = a->maxnz;
1135   c->solve_work       = 0;
1136   *B = C;
1137   return 0;
1138 }
1139 
1140 int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
1141 {
1142   Mat_SeqAIJ   *a;
1143   Mat          B;
1144   int          i, nz, ierr, fd, header[4],numtid,*rowlengths = 0,M,N,shift;
1145   PetscObject  vobj = (PetscObject) bview;
1146   MPI_Comm     comm = vobj->comm;
1147 
1148   MPI_Comm_size(comm,&numtid);
1149   if (numtid > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
1150   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1151   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
1152   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
1153   M = header[1]; N = header[2]; nz = header[3];
1154 
1155   /* read in row lengths */
1156   rowlengths = (int*) PETSCMALLOC( M*sizeof(int) ); CHKPTRQ(rowlengths);
1157   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
1158 
1159   /* create our matrix */
1160   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1161   B = *A;
1162   a = (Mat_SeqAIJ *) B->data;
1163   shift = a->indexshift;
1164 
1165   /* read in column indices and adjust for Fortran indexing*/
1166   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
1167   if (shift) {
1168     for ( i=0; i<nz; i++ ) {
1169       a->j[i] += 1;
1170     }
1171   }
1172 
1173   /* read in nonzero values */
1174   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
1175 
1176   /* set matrix "i" values */
1177   a->i[0] = -shift;
1178   for ( i=1; i<= M; i++ ) {
1179     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1180     a->ilen[i-1] = rowlengths[i-1];
1181   }
1182   PETSCFREE(rowlengths);
1183 
1184   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1185   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1186   return 0;
1187 }
1188 
1189 
1190 
1191