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