xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision 70c3da92039d8a9340517022d80c831ba45a21e4)
1*70c3da92SBarry Smith 
2*70c3da92SBarry Smith 
3*70c3da92SBarry Smith #ifndef lint
4*70c3da92SBarry Smith static char vcid[] = "$Id: aij.c,v 1.187 1996/09/23 16:21:45 bsmith Exp bsmith $";
5*70c3da92SBarry Smith #endif
6*70c3da92SBarry Smith 
7*70c3da92SBarry Smith /*
8*70c3da92SBarry Smith B    Defines the basic matrix operations for the AIJ (compressed row)
9*70c3da92SBarry Smith   matrix storage format.
10*70c3da92SBarry Smith */
11*70c3da92SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
12*70c3da92SBarry Smith #include "src/vec/vecimpl.h"
13*70c3da92SBarry Smith #include "src/inline/spops.h"
14*70c3da92SBarry Smith #include "petsc.h"
15*70c3da92SBarry Smith #include "src/inline/bitarray.h"
16*70c3da92SBarry Smith 
17*70c3da92SBarry Smith /*
18*70c3da92SBarry Smith     Basic AIJ format ILU based on drop tolerance
19*70c3da92SBarry Smith */
20*70c3da92SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
21*70c3da92SBarry Smith {
22*70c3da92SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
23*70c3da92SBarry Smith   int        ierr = 1;
24*70c3da92SBarry Smith 
25*70c3da92SBarry Smith   SETERRQ(ierr,"MatILUDTFactor_SeqAIJ:Not implemented");
26*70c3da92SBarry Smith }
27*70c3da92SBarry Smith 
28*70c3da92SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
29*70c3da92SBarry Smith 
30*70c3da92SBarry Smith static int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
31*70c3da92SBarry Smith                            PetscTruth *done)
32*70c3da92SBarry Smith {
33*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
34*70c3da92SBarry Smith   int        ierr,i,ishift;
35*70c3da92SBarry Smith 
36*70c3da92SBarry Smith   *n     = A->n;
37*70c3da92SBarry Smith   if (!ia) return 0;
38*70c3da92SBarry Smith   ishift = a->indexshift;
39*70c3da92SBarry Smith   if (symmetric) {
40*70c3da92SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
41*70c3da92SBarry Smith   } else if (oshift == 0 && ishift == -1) {
42*70c3da92SBarry Smith     int nz = a->i[a->n];
43*70c3da92SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
44*70c3da92SBarry Smith     *ia = (int *) PetscMalloc( (a->n+1)*sizeof(int) ); CHKPTRQ(*ia);
45*70c3da92SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
46*70c3da92SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
47*70c3da92SBarry Smith     for ( i=0; i<a->n+1; i++ ) (*ia)[i] = a->i[i] - 1;
48*70c3da92SBarry Smith   } else if (oshift == 1 && ishift == 0) {
49*70c3da92SBarry Smith     int nz = a->i[a->n] + 1;
50*70c3da92SBarry Smith     /* malloc space and  add 1 to i and j indices */
51*70c3da92SBarry Smith     *ia = (int *) PetscMalloc( (a->n+1)*sizeof(int) ); CHKPTRQ(*ia);
52*70c3da92SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
53*70c3da92SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
54*70c3da92SBarry Smith     for ( i=0; i<a->n+1; i++ ) (*ia)[i] = a->i[i] + 1;
55*70c3da92SBarry Smith   } else {
56*70c3da92SBarry Smith     *ia = a->i; *ja = a->j;
57*70c3da92SBarry Smith   }
58*70c3da92SBarry Smith 
59*70c3da92SBarry Smith   return 0;
60*70c3da92SBarry Smith }
61*70c3da92SBarry Smith 
62*70c3da92SBarry Smith static int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
63*70c3da92SBarry Smith                                PetscTruth *done)
64*70c3da92SBarry Smith {
65*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
66*70c3da92SBarry Smith   int        ishift = a->indexshift;
67*70c3da92SBarry Smith 
68*70c3da92SBarry Smith   if (!ia) return 0;
69*70c3da92SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
70*70c3da92SBarry Smith     PetscFree(*ia);
71*70c3da92SBarry Smith     PetscFree(*ja);
72*70c3da92SBarry Smith   }
73*70c3da92SBarry Smith   return 0;
74*70c3da92SBarry Smith }
75*70c3da92SBarry Smith 
76*70c3da92SBarry Smith static int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
77*70c3da92SBarry Smith                            PetscTruth *done)
78*70c3da92SBarry Smith {
79*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
80*70c3da92SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n;
81*70c3da92SBarry Smith   int        nz = a->i[n]+ishift,row,*jj,m,col;
82*70c3da92SBarry Smith 
83*70c3da92SBarry Smith   *nn     = A->n;
84*70c3da92SBarry Smith   if (!ia) return 0;
85*70c3da92SBarry Smith   if (symmetric) {
86*70c3da92SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
87*70c3da92SBarry Smith   } else {
88*70c3da92SBarry Smith     collengths = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(collengths);
89*70c3da92SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
90*70c3da92SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
91*70c3da92SBarry Smith     cja        = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cja);
92*70c3da92SBarry Smith     jj = a->j;
93*70c3da92SBarry Smith     for ( i=0; i<nz; i++ ) {
94*70c3da92SBarry Smith       collengths[jj[i] + ishift]++;
95*70c3da92SBarry Smith     }
96*70c3da92SBarry Smith     cia[0] = oshift;
97*70c3da92SBarry Smith     for ( i=0; i<n; i++) {
98*70c3da92SBarry Smith       cia[i+1] = cia[i] + collengths[i];
99*70c3da92SBarry Smith     }
100*70c3da92SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
101*70c3da92SBarry Smith     jj = a->j;
102*70c3da92SBarry Smith     for ( row=0; row<n; row++ ) {
103*70c3da92SBarry Smith       m = a->i[row+1] - a->i[row];
104*70c3da92SBarry Smith       for ( i=0; i<m; i++ ) {
105*70c3da92SBarry Smith 
106*70c3da92SBarry Smith         col = *jj++ + ishift;
107*70c3da92SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
108*70c3da92SBarry Smith       }
109*70c3da92SBarry Smith     }
110*70c3da92SBarry Smith     PetscFree(collengths);
111*70c3da92SBarry Smith     *ia = cia; *ja = cja;
112*70c3da92SBarry Smith   }
113*70c3da92SBarry Smith 
114*70c3da92SBarry Smith   return 0;
115*70c3da92SBarry Smith }
116*70c3da92SBarry Smith 
117*70c3da92SBarry Smith static int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
118*70c3da92SBarry Smith                                      int **ja,PetscTruth *done)
119*70c3da92SBarry Smith {
120*70c3da92SBarry Smith   if (!ia) return 0;
121*70c3da92SBarry Smith 
122*70c3da92SBarry Smith   PetscFree(*ia);
123*70c3da92SBarry Smith   PetscFree(*ja);
124*70c3da92SBarry Smith 
125*70c3da92SBarry Smith   return 0;
126*70c3da92SBarry Smith }
127*70c3da92SBarry Smith 
128*70c3da92SBarry Smith #define CHUNKSIZE   15
129*70c3da92SBarry Smith 
130*70c3da92SBarry Smith /* This version has row oriented v  */
131*70c3da92SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
132*70c3da92SBarry Smith {
133*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
134*70c3da92SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
135*70c3da92SBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
136*70c3da92SBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
137*70c3da92SBarry Smith   Scalar     *ap,value, *aa = a->a;
138*70c3da92SBarry Smith 
139*70c3da92SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
140*70c3da92SBarry Smith     row  = im[k];
141*70c3da92SBarry Smith #if defined(PETSC_BOPT_g)
142*70c3da92SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
143*70c3da92SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
144*70c3da92SBarry Smith #endif
145*70c3da92SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
146*70c3da92SBarry Smith     rmax = imax[row]; nrow = ailen[row];
147*70c3da92SBarry Smith     low = 0;
148*70c3da92SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
149*70c3da92SBarry Smith #if defined(PETSC_BOPT_g)
150*70c3da92SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
151*70c3da92SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
152*70c3da92SBarry Smith #endif
153*70c3da92SBarry Smith       col = in[l] - shift;
154*70c3da92SBarry Smith       if (roworiented) {
155*70c3da92SBarry Smith         value = *v++;
156*70c3da92SBarry Smith       }
157*70c3da92SBarry Smith       else {
158*70c3da92SBarry Smith         value = v[k + l*m];
159*70c3da92SBarry Smith       }
160*70c3da92SBarry Smith       if (!sorted) low = 0; high = nrow;
161*70c3da92SBarry Smith       while (high-low > 5) {
162*70c3da92SBarry Smith         t = (low+high)/2;
163*70c3da92SBarry Smith         if (rp[t] > col) high = t;
164*70c3da92SBarry Smith         else             low  = t;
165*70c3da92SBarry Smith       }
166*70c3da92SBarry Smith       for ( i=low; i<high; i++ ) {
167*70c3da92SBarry Smith         if (rp[i] > col) break;
168*70c3da92SBarry Smith         if (rp[i] == col) {
169*70c3da92SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
170*70c3da92SBarry Smith           else                  ap[i] = value;
171*70c3da92SBarry Smith           goto noinsert;
172*70c3da92SBarry Smith         }
173*70c3da92SBarry Smith       }
174*70c3da92SBarry Smith       if (nonew) goto noinsert;
175*70c3da92SBarry Smith       if (nrow >= rmax) {
176*70c3da92SBarry Smith         /* there is no extra room in row, therefore enlarge */
177*70c3da92SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
178*70c3da92SBarry Smith         Scalar *new_a;
179*70c3da92SBarry Smith 
180*70c3da92SBarry Smith         /* malloc new storage space */
181*70c3da92SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
182*70c3da92SBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
183*70c3da92SBarry Smith         new_j   = (int *) (new_a + new_nz);
184*70c3da92SBarry Smith         new_i   = new_j + new_nz;
185*70c3da92SBarry Smith 
186*70c3da92SBarry Smith         /* copy over old data into new slots */
187*70c3da92SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
188*70c3da92SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
189*70c3da92SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
190*70c3da92SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
191*70c3da92SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
192*70c3da92SBarry Smith                                                            len*sizeof(int));
193*70c3da92SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
194*70c3da92SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
195*70c3da92SBarry Smith                                                            len*sizeof(Scalar));
196*70c3da92SBarry Smith         /* free up old matrix storage */
197*70c3da92SBarry Smith         PetscFree(a->a);
198*70c3da92SBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
199*70c3da92SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
200*70c3da92SBarry Smith         a->singlemalloc = 1;
201*70c3da92SBarry Smith 
202*70c3da92SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
203*70c3da92SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
204*70c3da92SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
205*70c3da92SBarry Smith         a->maxnz += CHUNKSIZE;
206*70c3da92SBarry Smith         a->reallocs++;
207*70c3da92SBarry Smith       }
208*70c3da92SBarry Smith       N = nrow++ - 1; a->nz++;
209*70c3da92SBarry Smith       /* shift up all the later entries in this row */
210*70c3da92SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
211*70c3da92SBarry Smith         rp[ii+1] = rp[ii];
212*70c3da92SBarry Smith         ap[ii+1] = ap[ii];
213*70c3da92SBarry Smith       }
214*70c3da92SBarry Smith       rp[i] = col;
215*70c3da92SBarry Smith       ap[i] = value;
216*70c3da92SBarry Smith       noinsert:;
217*70c3da92SBarry Smith       low = i + 1;
218*70c3da92SBarry Smith     }
219*70c3da92SBarry Smith     ailen[row] = nrow;
220*70c3da92SBarry Smith   }
221*70c3da92SBarry Smith   return 0;
222*70c3da92SBarry Smith }
223*70c3da92SBarry Smith 
224*70c3da92SBarry Smith static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
225*70c3da92SBarry Smith {
226*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
227*70c3da92SBarry Smith   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
228*70c3da92SBarry Smith   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
229*70c3da92SBarry Smith   Scalar     *ap, *aa = a->a, zero = 0.0;
230*70c3da92SBarry Smith 
231*70c3da92SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over rows */
232*70c3da92SBarry Smith     row  = im[k];
233*70c3da92SBarry Smith     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
234*70c3da92SBarry Smith     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
235*70c3da92SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
236*70c3da92SBarry Smith     nrow = ailen[row];
237*70c3da92SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over columns */
238*70c3da92SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
239*70c3da92SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
240*70c3da92SBarry Smith       col = in[l] - shift;
241*70c3da92SBarry Smith       high = nrow; low = 0; /* assume unsorted */
242*70c3da92SBarry Smith       while (high-low > 5) {
243*70c3da92SBarry Smith         t = (low+high)/2;
244*70c3da92SBarry Smith         if (rp[t] > col) high = t;
245*70c3da92SBarry Smith         else             low  = t;
246*70c3da92SBarry Smith       }
247*70c3da92SBarry Smith       for ( i=low; i<high; i++ ) {
248*70c3da92SBarry Smith         if (rp[i] > col) break;
249*70c3da92SBarry Smith         if (rp[i] == col) {
250*70c3da92SBarry Smith           *v++ = ap[i];
251*70c3da92SBarry Smith           goto finished;
252*70c3da92SBarry Smith         }
253*70c3da92SBarry Smith       }
254*70c3da92SBarry Smith       *v++ = zero;
255*70c3da92SBarry Smith       finished:;
256*70c3da92SBarry Smith     }
257*70c3da92SBarry Smith   }
258*70c3da92SBarry Smith   return 0;
259*70c3da92SBarry Smith }
260*70c3da92SBarry Smith 
261*70c3da92SBarry Smith #include "draw.h"
262*70c3da92SBarry Smith #include "pinclude/pviewer.h"
263*70c3da92SBarry Smith #include "sys.h"
264*70c3da92SBarry Smith 
265*70c3da92SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
266*70c3da92SBarry Smith {
267*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
268*70c3da92SBarry Smith   int        i, fd, *col_lens, ierr;
269*70c3da92SBarry Smith 
270*70c3da92SBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
271*70c3da92SBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
272*70c3da92SBarry Smith   col_lens[0] = MAT_COOKIE;
273*70c3da92SBarry Smith   col_lens[1] = a->m;
274*70c3da92SBarry Smith   col_lens[2] = a->n;
275*70c3da92SBarry Smith   col_lens[3] = a->nz;
276*70c3da92SBarry Smith 
277*70c3da92SBarry Smith   /* store lengths of each row and write (including header) to file */
278*70c3da92SBarry Smith   for ( i=0; i<a->m; i++ ) {
279*70c3da92SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
280*70c3da92SBarry Smith   }
281*70c3da92SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
282*70c3da92SBarry Smith   PetscFree(col_lens);
283*70c3da92SBarry Smith 
284*70c3da92SBarry Smith   /* store column indices (zero start index) */
285*70c3da92SBarry Smith   if (a->indexshift) {
286*70c3da92SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
287*70c3da92SBarry Smith   }
288*70c3da92SBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr);
289*70c3da92SBarry Smith   if (a->indexshift) {
290*70c3da92SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
291*70c3da92SBarry Smith   }
292*70c3da92SBarry Smith 
293*70c3da92SBarry Smith   /* store nonzero values */
294*70c3da92SBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
295*70c3da92SBarry Smith   return 0;
296*70c3da92SBarry Smith }
297*70c3da92SBarry Smith 
298*70c3da92SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
299*70c3da92SBarry Smith {
300*70c3da92SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
301*70c3da92SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
302*70c3da92SBarry Smith   FILE        *fd;
303*70c3da92SBarry Smith   char        *outputname;
304*70c3da92SBarry Smith 
305*70c3da92SBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
306*70c3da92SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
307*70c3da92SBarry Smith   ierr = ViewerGetFormat(viewer,&format);
308*70c3da92SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
309*70c3da92SBarry Smith     return 0;
310*70c3da92SBarry Smith   }
311*70c3da92SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
312*70c3da92SBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
313*70c3da92SBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
314*70c3da92SBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
315*70c3da92SBarry Smith     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
316*70c3da92SBarry Smith         a->inode.node_count,a->inode.limit);
317*70c3da92SBarry Smith   }
318*70c3da92SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
319*70c3da92SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
320*70c3da92SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
321*70c3da92SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz);
322*70c3da92SBarry Smith     fprintf(fd,"zzz = [\n");
323*70c3da92SBarry Smith 
324*70c3da92SBarry Smith     for (i=0; i<m; i++) {
325*70c3da92SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
326*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
327*70c3da92SBarry Smith         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]),
328*70c3da92SBarry Smith                    imag(a->a[j]));
329*70c3da92SBarry Smith #else
330*70c3da92SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
331*70c3da92SBarry Smith #endif
332*70c3da92SBarry Smith       }
333*70c3da92SBarry Smith     }
334*70c3da92SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
335*70c3da92SBarry Smith   }
336*70c3da92SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
337*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
338*70c3da92SBarry Smith       fprintf(fd,"row %d:",i);
339*70c3da92SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
340*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
341*70c3da92SBarry Smith         if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0)
342*70c3da92SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
343*70c3da92SBarry Smith         else if (real(a->a[j]) != 0.0)
344*70c3da92SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
345*70c3da92SBarry Smith #else
346*70c3da92SBarry Smith         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
347*70c3da92SBarry Smith #endif
348*70c3da92SBarry Smith       }
349*70c3da92SBarry Smith       fprintf(fd,"\n");
350*70c3da92SBarry Smith     }
351*70c3da92SBarry Smith   }
352*70c3da92SBarry Smith   else {
353*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
354*70c3da92SBarry Smith       fprintf(fd,"row %d:",i);
355*70c3da92SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
356*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
357*70c3da92SBarry Smith         if (imag(a->a[j]) != 0.0) {
358*70c3da92SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
359*70c3da92SBarry Smith         }
360*70c3da92SBarry Smith         else {
361*70c3da92SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
362*70c3da92SBarry Smith         }
363*70c3da92SBarry Smith #else
364*70c3da92SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
365*70c3da92SBarry Smith #endif
366*70c3da92SBarry Smith       }
367*70c3da92SBarry Smith       fprintf(fd,"\n");
368*70c3da92SBarry Smith     }
369*70c3da92SBarry Smith   }
370*70c3da92SBarry Smith   fflush(fd);
371*70c3da92SBarry Smith   return 0;
372*70c3da92SBarry Smith }
373*70c3da92SBarry Smith 
374*70c3da92SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
375*70c3da92SBarry Smith {
376*70c3da92SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
377*70c3da92SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
378*70c3da92SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
379*70c3da92SBarry Smith   Draw        draw;
380*70c3da92SBarry Smith   DrawButton  button;
381*70c3da92SBarry Smith   PetscTruth  isnull;
382*70c3da92SBarry Smith 
383*70c3da92SBarry Smith   ViewerDrawGetDraw(viewer,&draw);
384*70c3da92SBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
385*70c3da92SBarry Smith 
386*70c3da92SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
387*70c3da92SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
388*70c3da92SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
389*70c3da92SBarry Smith   /* loop over matrix elements drawing boxes */
390*70c3da92SBarry Smith   color = DRAW_BLUE;
391*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
392*70c3da92SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
393*70c3da92SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
394*70c3da92SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
395*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
396*70c3da92SBarry Smith       if (real(a->a[j]) >=  0.) continue;
397*70c3da92SBarry Smith #else
398*70c3da92SBarry Smith       if (a->a[j] >=  0.) continue;
399*70c3da92SBarry Smith #endif
400*70c3da92SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
401*70c3da92SBarry Smith     }
402*70c3da92SBarry Smith   }
403*70c3da92SBarry Smith   color = DRAW_CYAN;
404*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
405*70c3da92SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
406*70c3da92SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
407*70c3da92SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
408*70c3da92SBarry Smith       if (a->a[j] !=  0.) continue;
409*70c3da92SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
410*70c3da92SBarry Smith     }
411*70c3da92SBarry Smith   }
412*70c3da92SBarry Smith   color = DRAW_RED;
413*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
414*70c3da92SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
415*70c3da92SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
416*70c3da92SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
417*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
418*70c3da92SBarry Smith       if (real(a->a[j]) <=  0.) continue;
419*70c3da92SBarry Smith #else
420*70c3da92SBarry Smith       if (a->a[j] <=  0.) continue;
421*70c3da92SBarry Smith #endif
422*70c3da92SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
423*70c3da92SBarry Smith     }
424*70c3da92SBarry Smith   }
425*70c3da92SBarry Smith   DrawFlush(draw);
426*70c3da92SBarry Smith   DrawGetPause(draw,&pause);
427*70c3da92SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
428*70c3da92SBarry Smith 
429*70c3da92SBarry Smith   /* allow the matrix to zoom or shrink */
430*70c3da92SBarry Smith   ierr = DrawCheckResizedWindow(draw);
431*70c3da92SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
432*70c3da92SBarry Smith   while (button != BUTTON_RIGHT) {
433*70c3da92SBarry Smith     DrawClear(draw);
434*70c3da92SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
435*70c3da92SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
436*70c3da92SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
437*70c3da92SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
438*70c3da92SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
439*70c3da92SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
440*70c3da92SBarry Smith     w *= scale; h *= scale;
441*70c3da92SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
442*70c3da92SBarry Smith     color = DRAW_BLUE;
443*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
444*70c3da92SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
445*70c3da92SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
446*70c3da92SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
447*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
448*70c3da92SBarry Smith         if (real(a->a[j]) >=  0.) continue;
449*70c3da92SBarry Smith #else
450*70c3da92SBarry Smith         if (a->a[j] >=  0.) continue;
451*70c3da92SBarry Smith #endif
452*70c3da92SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
453*70c3da92SBarry Smith       }
454*70c3da92SBarry Smith     }
455*70c3da92SBarry Smith     color = DRAW_CYAN;
456*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
457*70c3da92SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
458*70c3da92SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
459*70c3da92SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
460*70c3da92SBarry Smith         if (a->a[j] !=  0.) continue;
461*70c3da92SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
462*70c3da92SBarry Smith       }
463*70c3da92SBarry Smith     }
464*70c3da92SBarry Smith     color = DRAW_RED;
465*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
466*70c3da92SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
467*70c3da92SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
468*70c3da92SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
469*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
470*70c3da92SBarry Smith         if (real(a->a[j]) <=  0.) continue;
471*70c3da92SBarry Smith #else
472*70c3da92SBarry Smith         if (a->a[j] <=  0.) continue;
473*70c3da92SBarry Smith #endif
474*70c3da92SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
475*70c3da92SBarry Smith       }
476*70c3da92SBarry Smith     }
477*70c3da92SBarry Smith     ierr = DrawCheckResizedWindow(draw);
478*70c3da92SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
479*70c3da92SBarry Smith   }
480*70c3da92SBarry Smith   return 0;
481*70c3da92SBarry Smith }
482*70c3da92SBarry Smith 
483*70c3da92SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
484*70c3da92SBarry Smith {
485*70c3da92SBarry Smith   Mat         A = (Mat) obj;
486*70c3da92SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
487*70c3da92SBarry Smith   ViewerType  vtype;
488*70c3da92SBarry Smith   int         ierr;
489*70c3da92SBarry Smith 
490*70c3da92SBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
491*70c3da92SBarry Smith   if (vtype == MATLAB_VIEWER) {
492*70c3da92SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
493*70c3da92SBarry Smith   }
494*70c3da92SBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
495*70c3da92SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
496*70c3da92SBarry Smith   }
497*70c3da92SBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
498*70c3da92SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
499*70c3da92SBarry Smith   }
500*70c3da92SBarry Smith   else if (vtype == DRAW_VIEWER) {
501*70c3da92SBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
502*70c3da92SBarry Smith   }
503*70c3da92SBarry Smith   return 0;
504*70c3da92SBarry Smith }
505*70c3da92SBarry Smith 
506*70c3da92SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
507*70c3da92SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
508*70c3da92SBarry Smith {
509*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
510*70c3da92SBarry Smith   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
511*70c3da92SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
512*70c3da92SBarry Smith   Scalar     *aa = a->a, *ap;
513*70c3da92SBarry Smith 
514*70c3da92SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
515*70c3da92SBarry Smith 
516*70c3da92SBarry Smith   for ( i=1; i<m; i++ ) {
517*70c3da92SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
518*70c3da92SBarry Smith     fshift += imax[i-1] - ailen[i-1];
519*70c3da92SBarry Smith     if (fshift) {
520*70c3da92SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
521*70c3da92SBarry Smith       N = ailen[i];
522*70c3da92SBarry Smith       for ( j=0; j<N; j++ ) {
523*70c3da92SBarry Smith         ip[j-fshift] = ip[j];
524*70c3da92SBarry Smith         ap[j-fshift] = ap[j];
525*70c3da92SBarry Smith       }
526*70c3da92SBarry Smith     }
527*70c3da92SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
528*70c3da92SBarry Smith   }
529*70c3da92SBarry Smith   if (m) {
530*70c3da92SBarry Smith     fshift += imax[m-1] - ailen[m-1];
531*70c3da92SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
532*70c3da92SBarry Smith   }
533*70c3da92SBarry Smith   /* reset ilen and imax for each row */
534*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
535*70c3da92SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
536*70c3da92SBarry Smith   }
537*70c3da92SBarry Smith   a->nz = ai[m] + shift;
538*70c3da92SBarry Smith 
539*70c3da92SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
540*70c3da92SBarry Smith   if (fshift && a->diag) {
541*70c3da92SBarry Smith     PetscFree(a->diag);
542*70c3da92SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
543*70c3da92SBarry Smith     a->diag = 0;
544*70c3da92SBarry Smith   }
545*70c3da92SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
546*70c3da92SBarry Smith            m,a->n,fshift,a->nz);
547*70c3da92SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
548*70c3da92SBarry Smith            a->reallocs);
549*70c3da92SBarry Smith   A->info.nz_unneeded  = (double)fshift;
550*70c3da92SBarry Smith 
551*70c3da92SBarry Smith   /* check out for identical nodes. If found, use inode functions */
552*70c3da92SBarry Smith   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
553*70c3da92SBarry Smith   return 0;
554*70c3da92SBarry Smith }
555*70c3da92SBarry Smith 
556*70c3da92SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
557*70c3da92SBarry Smith {
558*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
559*70c3da92SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
560*70c3da92SBarry Smith   return 0;
561*70c3da92SBarry Smith }
562*70c3da92SBarry Smith 
563*70c3da92SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
564*70c3da92SBarry Smith {
565*70c3da92SBarry Smith   Mat        A  = (Mat) obj;
566*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
567*70c3da92SBarry Smith 
568*70c3da92SBarry Smith #if defined(PETSC_LOG)
569*70c3da92SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
570*70c3da92SBarry Smith #endif
571*70c3da92SBarry Smith   PetscFree(a->a);
572*70c3da92SBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
573*70c3da92SBarry Smith   if (a->diag) PetscFree(a->diag);
574*70c3da92SBarry Smith   if (a->ilen) PetscFree(a->ilen);
575*70c3da92SBarry Smith   if (a->imax) PetscFree(a->imax);
576*70c3da92SBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
577*70c3da92SBarry Smith   if (a->inode.size) PetscFree(a->inode.size);
578*70c3da92SBarry Smith   PetscFree(a);
579*70c3da92SBarry Smith   PLogObjectDestroy(A);
580*70c3da92SBarry Smith   PetscHeaderDestroy(A);
581*70c3da92SBarry Smith   return 0;
582*70c3da92SBarry Smith }
583*70c3da92SBarry Smith 
584*70c3da92SBarry Smith static int MatCompress_SeqAIJ(Mat A)
585*70c3da92SBarry Smith {
586*70c3da92SBarry Smith   return 0;
587*70c3da92SBarry Smith }
588*70c3da92SBarry Smith 
589*70c3da92SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
590*70c3da92SBarry Smith {
591*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
592*70c3da92SBarry Smith   if      (op == MAT_ROW_ORIENTED)              a->roworiented = 1;
593*70c3da92SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)           a->roworiented = 0;
594*70c3da92SBarry Smith   else if (op == MAT_COLUMNS_SORTED)            a->sorted      = 1;
595*70c3da92SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
596*70c3da92SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
597*70c3da92SBarry Smith   else if (op == MAT_ROWS_SORTED ||
598*70c3da92SBarry Smith            op == MAT_SYMMETRIC ||
599*70c3da92SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
600*70c3da92SBarry Smith            op == MAT_YES_NEW_DIAGONALS)
601*70c3da92SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
602*70c3da92SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
603*70c3da92SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:MAT_NO_NEW_DIAGONALS");}
604*70c3da92SBarry Smith   else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
605*70c3da92SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
606*70c3da92SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
607*70c3da92SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
608*70c3da92SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
609*70c3da92SBarry Smith   else
610*70c3da92SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
611*70c3da92SBarry Smith   return 0;
612*70c3da92SBarry Smith }
613*70c3da92SBarry Smith 
614*70c3da92SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
615*70c3da92SBarry Smith {
616*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
617*70c3da92SBarry Smith   int        i,j, n,shift = a->indexshift;
618*70c3da92SBarry Smith   Scalar     *x, zero = 0.0;
619*70c3da92SBarry Smith 
620*70c3da92SBarry Smith   VecSet(&zero,v);
621*70c3da92SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
622*70c3da92SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
623*70c3da92SBarry Smith   for ( i=0; i<a->m; i++ ) {
624*70c3da92SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
625*70c3da92SBarry Smith       if (a->j[j]+shift == i) {
626*70c3da92SBarry Smith         x[i] = a->a[j];
627*70c3da92SBarry Smith         break;
628*70c3da92SBarry Smith       }
629*70c3da92SBarry Smith     }
630*70c3da92SBarry Smith   }
631*70c3da92SBarry Smith   return 0;
632*70c3da92SBarry Smith }
633*70c3da92SBarry Smith 
634*70c3da92SBarry Smith /* -------------------------------------------------------*/
635*70c3da92SBarry Smith /* Should check that shapes of vectors and matrices match */
636*70c3da92SBarry Smith /* -------------------------------------------------------*/
637*70c3da92SBarry Smith int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
638*70c3da92SBarry Smith {
639*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
640*70c3da92SBarry Smith   Scalar     *x, *y, *v, alpha;
641*70c3da92SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
642*70c3da92SBarry Smith 
643*70c3da92SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
644*70c3da92SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
645*70c3da92SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
646*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
647*70c3da92SBarry Smith     idx   = a->j + a->i[i] + shift;
648*70c3da92SBarry Smith     v     = a->a + a->i[i] + shift;
649*70c3da92SBarry Smith     n     = a->i[i+1] - a->i[i];
650*70c3da92SBarry Smith     alpha = x[i];
651*70c3da92SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
652*70c3da92SBarry Smith   }
653*70c3da92SBarry Smith   PLogFlops(2*a->nz - a->n);
654*70c3da92SBarry Smith   return 0;
655*70c3da92SBarry Smith }
656*70c3da92SBarry Smith 
657*70c3da92SBarry Smith int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
658*70c3da92SBarry Smith {
659*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
660*70c3da92SBarry Smith   Scalar     *x, *y, *v, alpha;
661*70c3da92SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
662*70c3da92SBarry Smith 
663*70c3da92SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
664*70c3da92SBarry Smith   if (zz != yy) VecCopy(zz,yy);
665*70c3da92SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
666*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
667*70c3da92SBarry Smith     idx   = a->j + a->i[i] + shift;
668*70c3da92SBarry Smith     v     = a->a + a->i[i] + shift;
669*70c3da92SBarry Smith     n     = a->i[i+1] - a->i[i];
670*70c3da92SBarry Smith     alpha = x[i];
671*70c3da92SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
672*70c3da92SBarry Smith   }
673*70c3da92SBarry Smith   return 0;
674*70c3da92SBarry Smith }
675*70c3da92SBarry Smith 
676*70c3da92SBarry Smith int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
677*70c3da92SBarry Smith {
678*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
679*70c3da92SBarry Smith   Scalar     *x, *y, *v, sum;
680*70c3da92SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
681*70c3da92SBarry Smith 
682*70c3da92SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
683*70c3da92SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
684*70c3da92SBarry Smith   idx  = a->j;
685*70c3da92SBarry Smith   v    = a->a;
686*70c3da92SBarry Smith   ii   = a->i;
687*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
688*70c3da92SBarry Smith     n    = ii[1] - ii[0]; ii++;
689*70c3da92SBarry Smith     sum  = 0.0;
690*70c3da92SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
691*70c3da92SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
692*70c3da92SBarry Smith     while (n--) sum += *v++ * x[*idx++];
693*70c3da92SBarry Smith     y[i] = sum;
694*70c3da92SBarry Smith   }
695*70c3da92SBarry Smith   PLogFlops(2*a->nz - m);
696*70c3da92SBarry Smith   return 0;
697*70c3da92SBarry Smith }
698*70c3da92SBarry Smith 
699*70c3da92SBarry Smith int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
700*70c3da92SBarry Smith {
701*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
702*70c3da92SBarry Smith   Scalar     *x, *y, *z, *v, sum;
703*70c3da92SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
704*70c3da92SBarry Smith 
705*70c3da92SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
706*70c3da92SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
707*70c3da92SBarry Smith   idx  = a->j;
708*70c3da92SBarry Smith   v    = a->a;
709*70c3da92SBarry Smith   ii   = a->i;
710*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
711*70c3da92SBarry Smith     n    = ii[1] - ii[0]; ii++;
712*70c3da92SBarry Smith     sum  = y[i];
713*70c3da92SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
714*70c3da92SBarry Smith     while (n--) sum += *v++ * x[*idx++];
715*70c3da92SBarry Smith     z[i] = sum;
716*70c3da92SBarry Smith   }
717*70c3da92SBarry Smith   PLogFlops(2*a->nz);
718*70c3da92SBarry Smith   return 0;
719*70c3da92SBarry Smith }
720*70c3da92SBarry Smith 
721*70c3da92SBarry Smith /*
722*70c3da92SBarry Smith      Adds diagonal pointers to sparse matrix structure.
723*70c3da92SBarry Smith */
724*70c3da92SBarry Smith 
725*70c3da92SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
726*70c3da92SBarry Smith {
727*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
728*70c3da92SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
729*70c3da92SBarry Smith 
730*70c3da92SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
731*70c3da92SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
732*70c3da92SBarry Smith   for ( i=0; i<a->m; i++ ) {
733*70c3da92SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
734*70c3da92SBarry Smith       if (a->j[j]+shift == i) {
735*70c3da92SBarry Smith         diag[i] = j - shift;
736*70c3da92SBarry Smith         break;
737*70c3da92SBarry Smith       }
738*70c3da92SBarry Smith     }
739*70c3da92SBarry Smith   }
740*70c3da92SBarry Smith   a->diag = diag;
741*70c3da92SBarry Smith   return 0;
742*70c3da92SBarry Smith }
743*70c3da92SBarry Smith 
744*70c3da92SBarry Smith int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
745*70c3da92SBarry Smith                            double fshift,int its,Vec xx)
746*70c3da92SBarry Smith {
747*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
748*70c3da92SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
749*70c3da92SBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
750*70c3da92SBarry Smith 
751*70c3da92SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
752*70c3da92SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
753*70c3da92SBarry Smith   diag = a->diag;
754*70c3da92SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
755*70c3da92SBarry Smith   if (flag == SOR_APPLY_UPPER) {
756*70c3da92SBarry Smith    /* apply ( U + D/omega) to the vector */
757*70c3da92SBarry Smith     bs = b + shift;
758*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
759*70c3da92SBarry Smith         d    = fshift + a->a[diag[i] + shift];
760*70c3da92SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
761*70c3da92SBarry Smith         idx  = a->j + diag[i] + (!shift);
762*70c3da92SBarry Smith         v    = a->a + diag[i] + (!shift);
763*70c3da92SBarry Smith         sum  = b[i]*d/omega;
764*70c3da92SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
765*70c3da92SBarry Smith         x[i] = sum;
766*70c3da92SBarry Smith     }
767*70c3da92SBarry Smith     return 0;
768*70c3da92SBarry Smith   }
769*70c3da92SBarry Smith   if (flag == SOR_APPLY_LOWER) {
770*70c3da92SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
771*70c3da92SBarry Smith   }
772*70c3da92SBarry Smith   else if (flag & SOR_EISENSTAT) {
773*70c3da92SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
774*70c3da92SBarry Smith     U is upper triangular, E is diagonal; This routine applies
775*70c3da92SBarry Smith 
776*70c3da92SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
777*70c3da92SBarry Smith 
778*70c3da92SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
779*70c3da92SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
780*70c3da92SBarry Smith     is the relaxation factor.
781*70c3da92SBarry Smith     */
782*70c3da92SBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
783*70c3da92SBarry Smith     scale = (2.0/omega) - 1.0;
784*70c3da92SBarry Smith 
785*70c3da92SBarry Smith     /*  x = (E + U)^{-1} b */
786*70c3da92SBarry Smith     for ( i=m-1; i>=0; i-- ) {
787*70c3da92SBarry Smith       d    = fshift + a->a[diag[i] + shift];
788*70c3da92SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
789*70c3da92SBarry Smith       idx  = a->j + diag[i] + (!shift);
790*70c3da92SBarry Smith       v    = a->a + diag[i] + (!shift);
791*70c3da92SBarry Smith       sum  = b[i];
792*70c3da92SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
793*70c3da92SBarry Smith       x[i] = omega*(sum/d);
794*70c3da92SBarry Smith     }
795*70c3da92SBarry Smith 
796*70c3da92SBarry Smith     /*  t = b - (2*E - D)x */
797*70c3da92SBarry Smith     v = a->a;
798*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
799*70c3da92SBarry Smith 
800*70c3da92SBarry Smith     /*  t = (E + L)^{-1}t */
801*70c3da92SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
802*70c3da92SBarry Smith     diag = a->diag;
803*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
804*70c3da92SBarry Smith       d    = fshift + a->a[diag[i]+shift];
805*70c3da92SBarry Smith       n    = diag[i] - a->i[i];
806*70c3da92SBarry Smith       idx  = a->j + a->i[i] + shift;
807*70c3da92SBarry Smith       v    = a->a + a->i[i] + shift;
808*70c3da92SBarry Smith       sum  = t[i];
809*70c3da92SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
810*70c3da92SBarry Smith       t[i] = omega*(sum/d);
811*70c3da92SBarry Smith     }
812*70c3da92SBarry Smith 
813*70c3da92SBarry Smith     /*  x = x + t */
814*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
815*70c3da92SBarry Smith     PetscFree(t);
816*70c3da92SBarry Smith     return 0;
817*70c3da92SBarry Smith   }
818*70c3da92SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
819*70c3da92SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
820*70c3da92SBarry Smith       for ( i=0; i<m; i++ ) {
821*70c3da92SBarry Smith         d    = fshift + a->a[diag[i]+shift];
822*70c3da92SBarry Smith         n    = diag[i] - a->i[i];
823*70c3da92SBarry Smith         idx  = a->j + a->i[i] + shift;
824*70c3da92SBarry Smith         v    = a->a + a->i[i] + shift;
825*70c3da92SBarry Smith         sum  = b[i];
826*70c3da92SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
827*70c3da92SBarry Smith         x[i] = omega*(sum/d);
828*70c3da92SBarry Smith       }
829*70c3da92SBarry Smith       xb = x;
830*70c3da92SBarry Smith     }
831*70c3da92SBarry Smith     else xb = b;
832*70c3da92SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
833*70c3da92SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
834*70c3da92SBarry Smith       for ( i=0; i<m; i++ ) {
835*70c3da92SBarry Smith         x[i] *= a->a[diag[i]+shift];
836*70c3da92SBarry Smith       }
837*70c3da92SBarry Smith     }
838*70c3da92SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
839*70c3da92SBarry Smith       for ( i=m-1; i>=0; i-- ) {
840*70c3da92SBarry Smith         d    = fshift + a->a[diag[i] + shift];
841*70c3da92SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
842*70c3da92SBarry Smith         idx  = a->j + diag[i] + (!shift);
843*70c3da92SBarry Smith         v    = a->a + diag[i] + (!shift);
844*70c3da92SBarry Smith         sum  = xb[i];
845*70c3da92SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
846*70c3da92SBarry Smith         x[i] = omega*(sum/d);
847*70c3da92SBarry Smith       }
848*70c3da92SBarry Smith     }
849*70c3da92SBarry Smith     its--;
850*70c3da92SBarry Smith   }
851*70c3da92SBarry Smith   while (its--) {
852*70c3da92SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
853*70c3da92SBarry Smith       for ( i=0; i<m; i++ ) {
854*70c3da92SBarry Smith         d    = fshift + a->a[diag[i]+shift];
855*70c3da92SBarry Smith         n    = a->i[i+1] - a->i[i];
856*70c3da92SBarry Smith         idx  = a->j + a->i[i] + shift;
857*70c3da92SBarry Smith         v    = a->a + a->i[i] + shift;
858*70c3da92SBarry Smith         sum  = b[i];
859*70c3da92SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
860*70c3da92SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
861*70c3da92SBarry Smith       }
862*70c3da92SBarry Smith     }
863*70c3da92SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
864*70c3da92SBarry Smith       for ( i=m-1; i>=0; i-- ) {
865*70c3da92SBarry Smith         d    = fshift + a->a[diag[i] + shift];
866*70c3da92SBarry Smith         n    = a->i[i+1] - a->i[i];
867*70c3da92SBarry Smith         idx  = a->j + a->i[i] + shift;
868*70c3da92SBarry Smith         v    = a->a + a->i[i] + shift;
869*70c3da92SBarry Smith         sum  = b[i];
870*70c3da92SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
871*70c3da92SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
872*70c3da92SBarry Smith       }
873*70c3da92SBarry Smith     }
874*70c3da92SBarry Smith   }
875*70c3da92SBarry Smith   return 0;
876*70c3da92SBarry Smith }
877*70c3da92SBarry Smith 
878*70c3da92SBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
879*70c3da92SBarry Smith {
880*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
881*70c3da92SBarry Smith 
882*70c3da92SBarry Smith   info->rows_global    = (double)a->m;
883*70c3da92SBarry Smith   info->columns_global = (double)a->n;
884*70c3da92SBarry Smith   info->rows_local     = (double)a->m;
885*70c3da92SBarry Smith   info->columns_local  = (double)a->n;
886*70c3da92SBarry Smith   info->block_size     = 1.0;
887*70c3da92SBarry Smith   info->nz_allocated   = (double)a->maxnz;
888*70c3da92SBarry Smith   info->nz_used        = (double)a->nz;
889*70c3da92SBarry Smith   info->nz_unneeded    = (double)(a->maxnz - a->nz);
890*70c3da92SBarry Smith   /*  if (info->nz_unneeded != A->info.nz_unneeded)
891*70c3da92SBarry Smith     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
892*70c3da92SBarry Smith   info->assemblies     = (double)A->num_ass;
893*70c3da92SBarry Smith   info->mallocs        = (double)a->reallocs;
894*70c3da92SBarry Smith   info->memory         = A->mem;
895*70c3da92SBarry Smith   if (A->factor) {
896*70c3da92SBarry Smith     info->fill_ratio_given  = A->info.fill_ratio_given;
897*70c3da92SBarry Smith     info->fill_ratio_needed = A->info.fill_ratio_needed;
898*70c3da92SBarry Smith     info->factor_mallocs    = A->info.factor_mallocs;
899*70c3da92SBarry Smith   } else {
900*70c3da92SBarry Smith     info->fill_ratio_given  = 0;
901*70c3da92SBarry Smith     info->fill_ratio_needed = 0;
902*70c3da92SBarry Smith     info->factor_mallocs    = 0;
903*70c3da92SBarry Smith   }
904*70c3da92SBarry Smith   return 0;
905*70c3da92SBarry Smith }
906*70c3da92SBarry Smith 
907*70c3da92SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
908*70c3da92SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
909*70c3da92SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
910*70c3da92SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
911*70c3da92SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
912*70c3da92SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
913*70c3da92SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
914*70c3da92SBarry Smith 
915*70c3da92SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
916*70c3da92SBarry Smith {
917*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
918*70c3da92SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
919*70c3da92SBarry Smith 
920*70c3da92SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
921*70c3da92SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
922*70c3da92SBarry Smith   if (diag) {
923*70c3da92SBarry Smith     for ( i=0; i<N; i++ ) {
924*70c3da92SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
925*70c3da92SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
926*70c3da92SBarry Smith         a->ilen[rows[i]] = 1;
927*70c3da92SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
928*70c3da92SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
929*70c3da92SBarry Smith       }
930*70c3da92SBarry Smith       else {
931*70c3da92SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
932*70c3da92SBarry Smith         CHKERRQ(ierr);
933*70c3da92SBarry Smith       }
934*70c3da92SBarry Smith     }
935*70c3da92SBarry Smith   }
936*70c3da92SBarry Smith   else {
937*70c3da92SBarry Smith     for ( i=0; i<N; i++ ) {
938*70c3da92SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
939*70c3da92SBarry Smith       a->ilen[rows[i]] = 0;
940*70c3da92SBarry Smith     }
941*70c3da92SBarry Smith   }
942*70c3da92SBarry Smith   ISRestoreIndices(is,&rows);
943*70c3da92SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
944*70c3da92SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
945*70c3da92SBarry Smith   return 0;
946*70c3da92SBarry Smith }
947*70c3da92SBarry Smith 
948*70c3da92SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
949*70c3da92SBarry Smith {
950*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
951*70c3da92SBarry Smith   *m = a->m; *n = a->n;
952*70c3da92SBarry Smith   return 0;
953*70c3da92SBarry Smith }
954*70c3da92SBarry Smith 
955*70c3da92SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
956*70c3da92SBarry Smith {
957*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
958*70c3da92SBarry Smith   *m = 0; *n = a->m;
959*70c3da92SBarry Smith   return 0;
960*70c3da92SBarry Smith }
961*70c3da92SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
962*70c3da92SBarry Smith {
963*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
964*70c3da92SBarry Smith   int        *itmp,i,shift = a->indexshift;
965*70c3da92SBarry Smith 
966*70c3da92SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
967*70c3da92SBarry Smith 
968*70c3da92SBarry Smith   *nz = a->i[row+1] - a->i[row];
969*70c3da92SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
970*70c3da92SBarry Smith   if (idx) {
971*70c3da92SBarry Smith     itmp = a->j + a->i[row] + shift;
972*70c3da92SBarry Smith     if (*nz && shift) {
973*70c3da92SBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
974*70c3da92SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
975*70c3da92SBarry Smith     } else if (*nz) {
976*70c3da92SBarry Smith       *idx = itmp;
977*70c3da92SBarry Smith     }
978*70c3da92SBarry Smith     else *idx = 0;
979*70c3da92SBarry Smith   }
980*70c3da92SBarry Smith   return 0;
981*70c3da92SBarry Smith }
982*70c3da92SBarry Smith 
983*70c3da92SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
984*70c3da92SBarry Smith {
985*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
986*70c3da92SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
987*70c3da92SBarry Smith   return 0;
988*70c3da92SBarry Smith }
989*70c3da92SBarry Smith 
990*70c3da92SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
991*70c3da92SBarry Smith {
992*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
993*70c3da92SBarry Smith   Scalar     *v = a->a;
994*70c3da92SBarry Smith   double     sum = 0.0;
995*70c3da92SBarry Smith   int        i, j,shift = a->indexshift;
996*70c3da92SBarry Smith 
997*70c3da92SBarry Smith   if (type == NORM_FROBENIUS) {
998*70c3da92SBarry Smith     for (i=0; i<a->nz; i++ ) {
999*70c3da92SBarry Smith #if defined(PETSC_COMPLEX)
1000*70c3da92SBarry Smith       sum += real(conj(*v)*(*v)); v++;
1001*70c3da92SBarry Smith #else
1002*70c3da92SBarry Smith       sum += (*v)*(*v); v++;
1003*70c3da92SBarry Smith #endif
1004*70c3da92SBarry Smith     }
1005*70c3da92SBarry Smith     *norm = sqrt(sum);
1006*70c3da92SBarry Smith   }
1007*70c3da92SBarry Smith   else if (type == NORM_1) {
1008*70c3da92SBarry Smith     double *tmp;
1009*70c3da92SBarry Smith     int    *jj = a->j;
1010*70c3da92SBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
1011*70c3da92SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
1012*70c3da92SBarry Smith     *norm = 0.0;
1013*70c3da92SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1014*70c3da92SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
1015*70c3da92SBarry Smith     }
1016*70c3da92SBarry Smith     for ( j=0; j<a->n; j++ ) {
1017*70c3da92SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
1018*70c3da92SBarry Smith     }
1019*70c3da92SBarry Smith     PetscFree(tmp);
1020*70c3da92SBarry Smith   }
1021*70c3da92SBarry Smith   else if (type == NORM_INFINITY) {
1022*70c3da92SBarry Smith     *norm = 0.0;
1023*70c3da92SBarry Smith     for ( j=0; j<a->m; j++ ) {
1024*70c3da92SBarry Smith       v = a->a + a->i[j] + shift;
1025*70c3da92SBarry Smith       sum = 0.0;
1026*70c3da92SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1027*70c3da92SBarry Smith         sum += PetscAbsScalar(*v); v++;
1028*70c3da92SBarry Smith       }
1029*70c3da92SBarry Smith       if (sum > *norm) *norm = sum;
1030*70c3da92SBarry Smith     }
1031*70c3da92SBarry Smith   }
1032*70c3da92SBarry Smith   else {
1033*70c3da92SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
1034*70c3da92SBarry Smith   }
1035*70c3da92SBarry Smith   return 0;
1036*70c3da92SBarry Smith }
1037*70c3da92SBarry Smith 
1038*70c3da92SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
1039*70c3da92SBarry Smith {
1040*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1041*70c3da92SBarry Smith   Mat        C;
1042*70c3da92SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1043*70c3da92SBarry Smith   int        shift = a->indexshift;
1044*70c3da92SBarry Smith   Scalar     *array = a->a;
1045*70c3da92SBarry Smith 
1046*70c3da92SBarry Smith   if (B == PETSC_NULL && m != a->n)
1047*70c3da92SBarry Smith     SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place");
1048*70c3da92SBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1049*70c3da92SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
1050*70c3da92SBarry Smith   if (shift) {
1051*70c3da92SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
1052*70c3da92SBarry Smith   }
1053*70c3da92SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1054*70c3da92SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
1055*70c3da92SBarry Smith   PetscFree(col);
1056*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
1057*70c3da92SBarry Smith     len = ai[i+1]-ai[i];
1058*70c3da92SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
1059*70c3da92SBarry Smith     array += len; aj += len;
1060*70c3da92SBarry Smith   }
1061*70c3da92SBarry Smith   if (shift) {
1062*70c3da92SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
1063*70c3da92SBarry Smith   }
1064*70c3da92SBarry Smith 
1065*70c3da92SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1066*70c3da92SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1067*70c3da92SBarry Smith 
1068*70c3da92SBarry Smith   if (B != PETSC_NULL) {
1069*70c3da92SBarry Smith     *B = C;
1070*70c3da92SBarry Smith   } else {
1071*70c3da92SBarry Smith     /* This isn't really an in-place transpose */
1072*70c3da92SBarry Smith     PetscFree(a->a);
1073*70c3da92SBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
1074*70c3da92SBarry Smith     if (a->diag) PetscFree(a->diag);
1075*70c3da92SBarry Smith     if (a->ilen) PetscFree(a->ilen);
1076*70c3da92SBarry Smith     if (a->imax) PetscFree(a->imax);
1077*70c3da92SBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
1078*70c3da92SBarry Smith     if (a->inode.size) PetscFree(a->inode.size);
1079*70c3da92SBarry Smith     PetscFree(a);
1080*70c3da92SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
1081*70c3da92SBarry Smith     PetscHeaderDestroy(C);
1082*70c3da92SBarry Smith   }
1083*70c3da92SBarry Smith   return 0;
1084*70c3da92SBarry Smith }
1085*70c3da92SBarry Smith 
1086*70c3da92SBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1087*70c3da92SBarry Smith {
1088*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1089*70c3da92SBarry Smith   Scalar     *l,*r,x,*v;
1090*70c3da92SBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
1091*70c3da92SBarry Smith 
1092*70c3da92SBarry Smith   if (ll) {
1093*70c3da92SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
1094*70c3da92SBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
1095*70c3da92SBarry Smith     v = a->a;
1096*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
1097*70c3da92SBarry Smith       x = l[i];
1098*70c3da92SBarry Smith       M = a->i[i+1] - a->i[i];
1099*70c3da92SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
1100*70c3da92SBarry Smith     }
1101*70c3da92SBarry Smith     PLogFlops(nz);
1102*70c3da92SBarry Smith   }
1103*70c3da92SBarry Smith   if (rr) {
1104*70c3da92SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
1105*70c3da92SBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
1106*70c3da92SBarry Smith     v = a->a; jj = a->j;
1107*70c3da92SBarry Smith     for ( i=0; i<nz; i++ ) {
1108*70c3da92SBarry Smith       (*v++) *= r[*jj++ + shift];
1109*70c3da92SBarry Smith     }
1110*70c3da92SBarry Smith     PLogFlops(nz);
1111*70c3da92SBarry Smith   }
1112*70c3da92SBarry Smith   return 0;
1113*70c3da92SBarry Smith }
1114*70c3da92SBarry Smith 
1115*70c3da92SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
1116*70c3da92SBarry Smith {
1117*70c3da92SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1118*70c3da92SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
1119*70c3da92SBarry Smith   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1120*70c3da92SBarry Smith   register int sum,lensi;
1121*70c3da92SBarry Smith   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
1122*70c3da92SBarry Smith   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
1123*70c3da92SBarry Smith   Scalar       *a_new,*mat_a;
1124*70c3da92SBarry Smith   Mat          C;
1125*70c3da92SBarry Smith 
1126*70c3da92SBarry Smith   ierr = ISSorted(isrow,(PetscTruth*)&i);
1127*70c3da92SBarry Smith   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:ISrow is not sorted");
1128*70c3da92SBarry Smith   ierr = ISSorted(iscol,(PetscTruth*)&i);
1129*70c3da92SBarry Smith   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IScol is not sorted");
1130*70c3da92SBarry Smith 
1131*70c3da92SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
1132*70c3da92SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
1133*70c3da92SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
1134*70c3da92SBarry Smith 
1135*70c3da92SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
1136*70c3da92SBarry Smith     /* special case of contiguous rows */
1137*70c3da92SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
1138*70c3da92SBarry Smith     starts = lens + ncols;
1139*70c3da92SBarry Smith     /* loop over new rows determining lens and starting points */
1140*70c3da92SBarry Smith     for (i=0; i<nrows; i++) {
1141*70c3da92SBarry Smith       kstart  = ai[irow[i]]+shift;
1142*70c3da92SBarry Smith       kend    = kstart + ailen[irow[i]];
1143*70c3da92SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1144*70c3da92SBarry Smith         if (aj[k]+shift >= first) {
1145*70c3da92SBarry Smith           starts[i] = k;
1146*70c3da92SBarry Smith           break;
1147*70c3da92SBarry Smith 	}
1148*70c3da92SBarry Smith       }
1149*70c3da92SBarry Smith       sum = 0;
1150*70c3da92SBarry Smith       while (k < kend) {
1151*70c3da92SBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1152*70c3da92SBarry Smith         sum++;
1153*70c3da92SBarry Smith       }
1154*70c3da92SBarry Smith       lens[i] = sum;
1155*70c3da92SBarry Smith     }
1156*70c3da92SBarry Smith     /* create submatrix */
1157*70c3da92SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
1158*70c3da92SBarry Smith       int n_cols,n_rows;
1159*70c3da92SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1160*70c3da92SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1161*70c3da92SBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
1162*70c3da92SBarry Smith       C = *B;
1163*70c3da92SBarry Smith     }
1164*70c3da92SBarry Smith     else {
1165*70c3da92SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
1166*70c3da92SBarry Smith     }
1167*70c3da92SBarry Smith     c = (Mat_SeqAIJ*) C->data;
1168*70c3da92SBarry Smith 
1169*70c3da92SBarry Smith     /* loop over rows inserting into submatrix */
1170*70c3da92SBarry Smith     a_new    = c->a;
1171*70c3da92SBarry Smith     j_new    = c->j;
1172*70c3da92SBarry Smith     i_new    = c->i;
1173*70c3da92SBarry Smith     i_new[0] = -shift;
1174*70c3da92SBarry Smith     for (i=0; i<nrows; i++) {
1175*70c3da92SBarry Smith       ii    = starts[i];
1176*70c3da92SBarry Smith       lensi = lens[i];
1177*70c3da92SBarry Smith       for ( k=0; k<lensi; k++ ) {
1178*70c3da92SBarry Smith         *j_new++ = aj[ii+k] - first;
1179*70c3da92SBarry Smith       }
1180*70c3da92SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1181*70c3da92SBarry Smith       a_new      += lensi;
1182*70c3da92SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1183*70c3da92SBarry Smith       c->ilen[i]  = lensi;
1184*70c3da92SBarry Smith     }
1185*70c3da92SBarry Smith     PetscFree(lens);
1186*70c3da92SBarry Smith   }
1187*70c3da92SBarry Smith   else {
1188*70c3da92SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
1189*70c3da92SBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
1190*70c3da92SBarry Smith     ssmap = smap + shift;
1191*70c3da92SBarry Smith     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1192*70c3da92SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
1193*70c3da92SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
1194*70c3da92SBarry Smith     /* determine lens of each row */
1195*70c3da92SBarry Smith     for (i=0; i<nrows; i++) {
1196*70c3da92SBarry Smith       kstart  = ai[irow[i]]+shift;
1197*70c3da92SBarry Smith       kend    = kstart + a->ilen[irow[i]];
1198*70c3da92SBarry Smith       lens[i] = 0;
1199*70c3da92SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1200*70c3da92SBarry Smith         if (ssmap[aj[k]]) {
1201*70c3da92SBarry Smith           lens[i]++;
1202*70c3da92SBarry Smith         }
1203*70c3da92SBarry Smith       }
1204*70c3da92SBarry Smith     }
1205*70c3da92SBarry Smith     /* Create and fill new matrix */
1206*70c3da92SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
1207*70c3da92SBarry Smith       c = (Mat_SeqAIJ *)((*B)->data);
1208*70c3da92SBarry Smith 
1209*70c3da92SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1210*70c3da92SBarry Smith       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1211*70c3da92SBarry Smith         SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros");
1212*70c3da92SBarry Smith       }
1213*70c3da92SBarry Smith       PetscMemzero(c->ilen,c->m*sizeof(int));
1214*70c3da92SBarry Smith       C = *B;
1215*70c3da92SBarry Smith     }
1216*70c3da92SBarry Smith     else {
1217*70c3da92SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
1218*70c3da92SBarry Smith     }
1219*70c3da92SBarry Smith     c = (Mat_SeqAIJ *)(C->data);
1220*70c3da92SBarry Smith     for (i=0; i<nrows; i++) {
1221*70c3da92SBarry Smith       row    = irow[i];
1222*70c3da92SBarry Smith       nznew  = 0;
1223*70c3da92SBarry Smith       kstart = ai[row]+shift;
1224*70c3da92SBarry Smith       kend   = kstart + a->ilen[row];
1225*70c3da92SBarry Smith       mat_i  = c->i[i]+shift;
1226*70c3da92SBarry Smith       mat_j  = c->j + mat_i;
1227*70c3da92SBarry Smith       mat_a  = c->a + mat_i;
1228*70c3da92SBarry Smith       mat_ilen = c->ilen + i;
1229*70c3da92SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1230*70c3da92SBarry Smith         if ((tcol=ssmap[a->j[k]])) {
1231*70c3da92SBarry Smith           *mat_j++ = tcol - (!shift);
1232*70c3da92SBarry Smith           *mat_a++ = a->a[k];
1233*70c3da92SBarry Smith           (*mat_ilen)++;
1234*70c3da92SBarry Smith 
1235*70c3da92SBarry Smith         }
1236*70c3da92SBarry Smith       }
1237*70c3da92SBarry Smith     }
1238*70c3da92SBarry Smith     /* Free work space */
1239*70c3da92SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1240*70c3da92SBarry Smith     PetscFree(smap); PetscFree(lens);
1241*70c3da92SBarry Smith   }
1242*70c3da92SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1243*70c3da92SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1244*70c3da92SBarry Smith 
1245*70c3da92SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1246*70c3da92SBarry Smith   *B = C;
1247*70c3da92SBarry Smith   return 0;
1248*70c3da92SBarry Smith }
1249*70c3da92SBarry Smith 
1250*70c3da92SBarry Smith /*
1251*70c3da92SBarry Smith      note: This can only work for identity for row and col. It would
1252*70c3da92SBarry Smith    be good to check this and otherwise generate an error.
1253*70c3da92SBarry Smith */
1254*70c3da92SBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1255*70c3da92SBarry Smith {
1256*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1257*70c3da92SBarry Smith   int        ierr;
1258*70c3da92SBarry Smith   Mat        outA;
1259*70c3da92SBarry Smith 
1260*70c3da92SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1261*70c3da92SBarry Smith 
1262*70c3da92SBarry Smith   outA          = inA;
1263*70c3da92SBarry Smith   inA->factor   = FACTOR_LU;
1264*70c3da92SBarry Smith   a->row        = row;
1265*70c3da92SBarry Smith   a->col        = col;
1266*70c3da92SBarry Smith 
1267*70c3da92SBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
1268*70c3da92SBarry Smith 
1269*70c3da92SBarry Smith   if (!a->diag) {
1270*70c3da92SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
1271*70c3da92SBarry Smith   }
1272*70c3da92SBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1273*70c3da92SBarry Smith   return 0;
1274*70c3da92SBarry Smith }
1275*70c3da92SBarry Smith 
1276*70c3da92SBarry Smith #include "pinclude/plapack.h"
1277*70c3da92SBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1278*70c3da92SBarry Smith {
1279*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1280*70c3da92SBarry Smith   int        one = 1;
1281*70c3da92SBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1282*70c3da92SBarry Smith   PLogFlops(a->nz);
1283*70c3da92SBarry Smith   return 0;
1284*70c3da92SBarry Smith }
1285*70c3da92SBarry Smith 
1286*70c3da92SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1287*70c3da92SBarry Smith                                     Mat **B)
1288*70c3da92SBarry Smith {
1289*70c3da92SBarry Smith   int ierr,i;
1290*70c3da92SBarry Smith 
1291*70c3da92SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1292*70c3da92SBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1293*70c3da92SBarry Smith   }
1294*70c3da92SBarry Smith 
1295*70c3da92SBarry Smith   for ( i=0; i<n; i++ ) {
1296*70c3da92SBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1297*70c3da92SBarry Smith   }
1298*70c3da92SBarry Smith   return 0;
1299*70c3da92SBarry Smith }
1300*70c3da92SBarry Smith 
1301*70c3da92SBarry Smith static int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
1302*70c3da92SBarry Smith {
1303*70c3da92SBarry Smith   *bs = 1;
1304*70c3da92SBarry Smith   return 0;
1305*70c3da92SBarry Smith }
1306*70c3da92SBarry Smith 
1307*70c3da92SBarry Smith static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
1308*70c3da92SBarry Smith {
1309*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1310*70c3da92SBarry Smith   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
1311*70c3da92SBarry Smith   int        start, end, *ai, *aj;
1312*70c3da92SBarry Smith   char       *table;
1313*70c3da92SBarry Smith   shift = a->indexshift;
1314*70c3da92SBarry Smith   m     = a->m;
1315*70c3da92SBarry Smith   ai    = a->i;
1316*70c3da92SBarry Smith   aj    = a->j+shift;
1317*70c3da92SBarry Smith 
1318*70c3da92SBarry Smith   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
1319*70c3da92SBarry Smith 
1320*70c3da92SBarry Smith   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
1321*70c3da92SBarry Smith   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1322*70c3da92SBarry Smith 
1323*70c3da92SBarry Smith   for ( i=0; i<is_max; i++ ) {
1324*70c3da92SBarry Smith     /* Initialise the two local arrays */
1325*70c3da92SBarry Smith     isz  = 0;
1326*70c3da92SBarry Smith     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1327*70c3da92SBarry Smith 
1328*70c3da92SBarry Smith                 /* Extract the indices, assume there can be duplicate entries */
1329*70c3da92SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
1330*70c3da92SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1331*70c3da92SBarry Smith 
1332*70c3da92SBarry Smith     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1333*70c3da92SBarry Smith     for ( j=0; j<n ; ++j){
1334*70c3da92SBarry Smith       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
1335*70c3da92SBarry Smith     }
1336*70c3da92SBarry Smith     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
1337*70c3da92SBarry Smith     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1338*70c3da92SBarry Smith 
1339*70c3da92SBarry Smith     k = 0;
1340*70c3da92SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
1341*70c3da92SBarry Smith       n = isz;
1342*70c3da92SBarry Smith       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1343*70c3da92SBarry Smith         row   = nidx[k];
1344*70c3da92SBarry Smith         start = ai[row];
1345*70c3da92SBarry Smith         end   = ai[row+1];
1346*70c3da92SBarry Smith         for ( l = start; l<end ; l++){
1347*70c3da92SBarry Smith           val = aj[l] + shift;
1348*70c3da92SBarry Smith           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1349*70c3da92SBarry Smith         }
1350*70c3da92SBarry Smith       }
1351*70c3da92SBarry Smith     }
1352*70c3da92SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1353*70c3da92SBarry Smith   }
1354*70c3da92SBarry Smith   PetscFree(table);
1355*70c3da92SBarry Smith   PetscFree(nidx);
1356*70c3da92SBarry Smith   return 0;
1357*70c3da92SBarry Smith }
1358*70c3da92SBarry Smith 
1359*70c3da92SBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1360*70c3da92SBarry Smith {
1361*70c3da92SBarry Smith   static int called = 0;
1362*70c3da92SBarry Smith   MPI_Comm   comm = A->comm;
1363*70c3da92SBarry Smith 
1364*70c3da92SBarry Smith   if (called) return 0; else called = 1;
1365*70c3da92SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1366*70c3da92SBarry Smith   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
1367*70c3da92SBarry Smith   PetscPrintf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
1368*70c3da92SBarry Smith   PetscPrintf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
1369*70c3da92SBarry Smith   PetscPrintf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1370*70c3da92SBarry Smith #if defined(HAVE_ESSL)
1371*70c3da92SBarry Smith   PetscPrintf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1372*70c3da92SBarry Smith #endif
1373*70c3da92SBarry Smith   return 0;
1374*70c3da92SBarry Smith }
1375*70c3da92SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1376*70c3da92SBarry Smith /* -------------------------------------------------------------------*/
1377*70c3da92SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
1378*70c3da92SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1379*70c3da92SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1380*70c3da92SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
1381*70c3da92SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
1382*70c3da92SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
1383*70c3da92SBarry Smith        MatLUFactor_SeqAIJ,0,
1384*70c3da92SBarry Smith        MatRelax_SeqAIJ,
1385*70c3da92SBarry Smith        MatTranspose_SeqAIJ,
1386*70c3da92SBarry Smith        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1387*70c3da92SBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
1388*70c3da92SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
1389*70c3da92SBarry Smith        MatCompress_SeqAIJ,
1390*70c3da92SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
1391*70c3da92SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
1392*70c3da92SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
1393*70c3da92SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
1394*70c3da92SBarry Smith        0,0,MatConvert_SeqAIJ,
1395*70c3da92SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1396*70c3da92SBarry Smith        MatILUFactor_SeqAIJ,0,0,
1397*70c3da92SBarry Smith        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1398*70c3da92SBarry Smith        MatGetValues_SeqAIJ,0,
1399*70c3da92SBarry Smith        MatPrintHelp_SeqAIJ,
1400*70c3da92SBarry Smith        MatScale_SeqAIJ,0,0,
1401*70c3da92SBarry Smith        MatILUDTFactor_SeqAIJ,
1402*70c3da92SBarry Smith        MatGetBlockSize_SeqAIJ,
1403*70c3da92SBarry Smith        MatGetRowIJ_SeqAIJ,
1404*70c3da92SBarry Smith        MatRestoreRowIJ_SeqAIJ,
1405*70c3da92SBarry Smith        MatGetColumnIJ_SeqAIJ,
1406*70c3da92SBarry Smith        MatRestoreColumnIJ_SeqAIJ};
1407*70c3da92SBarry Smith 
1408*70c3da92SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
1409*70c3da92SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
1410*70c3da92SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
1411*70c3da92SBarry Smith 
1412*70c3da92SBarry Smith /*@C
1413*70c3da92SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
1414*70c3da92SBarry Smith    (the default parallel PETSc format).  For good matrix assembly performance
1415*70c3da92SBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
1416*70c3da92SBarry Smith    (or the array nzz).  By setting these parameters accurately, performance
1417*70c3da92SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
1418*70c3da92SBarry Smith 
1419*70c3da92SBarry Smith    Input Parameters:
1420*70c3da92SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
1421*70c3da92SBarry Smith .  m - number of rows
1422*70c3da92SBarry Smith .  n - number of columns
1423*70c3da92SBarry Smith .  nz - number of nonzeros per row (same for all rows)
1424*70c3da92SBarry Smith .  nzz - array containing the number of nonzeros in the various rows
1425*70c3da92SBarry Smith          (possibly different for each row) or PETSC_NULL
1426*70c3da92SBarry Smith 
1427*70c3da92SBarry Smith    Output Parameter:
1428*70c3da92SBarry Smith .  A - the matrix
1429*70c3da92SBarry Smith 
1430*70c3da92SBarry Smith    Notes:
1431*70c3da92SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
1432*70c3da92SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
1433*70c3da92SBarry Smith    storage.  That is, the stored row and column indices can begin at
1434*70c3da92SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1435*70c3da92SBarry Smith 
1436*70c3da92SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1437*70c3da92SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1438*70c3da92SBarry Smith    allocation.  For large problems you MUST preallocate memory or you
1439*70c3da92SBarry Smith    will get TERRIBLE performance, see the users' manual chapter on
1440*70c3da92SBarry Smith    matrices and the file $(PETSC_DIR)/Performance.
1441*70c3da92SBarry Smith 
1442*70c3da92SBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1443*70c3da92SBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1444*70c3da92SBarry Smith    search for consecutive rows with the same nonzero structure, thereby
1445*70c3da92SBarry Smith    reusing matrix information to achieve increased efficiency.
1446*70c3da92SBarry Smith 
1447*70c3da92SBarry Smith    Options Database Keys:
1448*70c3da92SBarry Smith $    -mat_aij_no_inode  - Do not use inodes
1449*70c3da92SBarry Smith $    -mat_aij_inode_limit <limit> - Set inode limit.
1450*70c3da92SBarry Smith $        (max limit=5)
1451*70c3da92SBarry Smith $    -mat_aij_oneindex - Internally use indexing starting at 1
1452*70c3da92SBarry Smith $        rather than 0.  Note: When calling MatSetValues(),
1453*70c3da92SBarry Smith $        the user still MUST index entries starting at 0!
1454*70c3da92SBarry Smith 
1455*70c3da92SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
1456*70c3da92SBarry Smith @*/
1457*70c3da92SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
1458*70c3da92SBarry Smith {
1459*70c3da92SBarry Smith   Mat        B;
1460*70c3da92SBarry Smith   Mat_SeqAIJ *b;
1461*70c3da92SBarry Smith   int        i, len, ierr, flg,size;
1462*70c3da92SBarry Smith 
1463*70c3da92SBarry Smith   MPI_Comm_size(comm,&size);
1464*70c3da92SBarry Smith   if (size > 1) SETERRQ(1,"MatCreateSeqAIJ:Comm must be of size 1");
1465*70c3da92SBarry Smith 
1466*70c3da92SBarry Smith   *A                  = 0;
1467*70c3da92SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1468*70c3da92SBarry Smith   PLogObjectCreate(B);
1469*70c3da92SBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1470*70c3da92SBarry Smith   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1471*70c3da92SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1472*70c3da92SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1473*70c3da92SBarry Smith   B->view             = MatView_SeqAIJ;
1474*70c3da92SBarry Smith   B->factor           = 0;
1475*70c3da92SBarry Smith   B->lupivotthreshold = 1.0;
1476*70c3da92SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
1477*70c3da92SBarry Smith                           &flg); CHKERRQ(ierr);
1478*70c3da92SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
1479*70c3da92SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
1480*70c3da92SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1481*70c3da92SBarry Smith   b->row              = 0;
1482*70c3da92SBarry Smith   b->col              = 0;
1483*70c3da92SBarry Smith   b->indexshift       = 0;
1484*70c3da92SBarry Smith   b->reallocs         = 0;
1485*70c3da92SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
1486*70c3da92SBarry Smith   if (flg) b->indexshift = -1;
1487*70c3da92SBarry Smith 
1488*70c3da92SBarry Smith   b->m = m; B->m = m; B->M = m;
1489*70c3da92SBarry Smith   b->n = n; B->n = n; B->N = n;
1490*70c3da92SBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1491*70c3da92SBarry Smith   if (nnz == PETSC_NULL) {
1492*70c3da92SBarry Smith     if (nz == PETSC_DEFAULT) nz = 10;
1493*70c3da92SBarry Smith     else if (nz <= 0)        nz = 1;
1494*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
1495*70c3da92SBarry Smith     nz = nz*m;
1496*70c3da92SBarry Smith   }
1497*70c3da92SBarry Smith   else {
1498*70c3da92SBarry Smith     nz = 0;
1499*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
1500*70c3da92SBarry Smith   }
1501*70c3da92SBarry Smith 
1502*70c3da92SBarry Smith   /* allocate the matrix space */
1503*70c3da92SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
1504*70c3da92SBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1505*70c3da92SBarry Smith   b->j  = (int *) (b->a + nz);
1506*70c3da92SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1507*70c3da92SBarry Smith   b->i  = b->j + nz;
1508*70c3da92SBarry Smith   b->singlemalloc = 1;
1509*70c3da92SBarry Smith 
1510*70c3da92SBarry Smith   b->i[0] = -b->indexshift;
1511*70c3da92SBarry Smith   for (i=1; i<m+1; i++) {
1512*70c3da92SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
1513*70c3da92SBarry Smith   }
1514*70c3da92SBarry Smith 
1515*70c3da92SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
1516*70c3da92SBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1517*70c3da92SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1518*70c3da92SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
1519*70c3da92SBarry Smith 
1520*70c3da92SBarry Smith   b->nz               = 0;
1521*70c3da92SBarry Smith   b->maxnz            = nz;
1522*70c3da92SBarry Smith   b->sorted           = 0;
1523*70c3da92SBarry Smith   b->roworiented      = 1;
1524*70c3da92SBarry Smith   b->nonew            = 0;
1525*70c3da92SBarry Smith   b->diag             = 0;
1526*70c3da92SBarry Smith   b->solve_work       = 0;
1527*70c3da92SBarry Smith   b->spptr            = 0;
1528*70c3da92SBarry Smith   b->inode.node_count = 0;
1529*70c3da92SBarry Smith   b->inode.size       = 0;
1530*70c3da92SBarry Smith   b->inode.limit      = 5;
1531*70c3da92SBarry Smith   b->inode.max_limit  = 5;
1532*70c3da92SBarry Smith   B->info.nz_unneeded = (double)b->maxnz;
1533*70c3da92SBarry Smith 
1534*70c3da92SBarry Smith   *A = B;
1535*70c3da92SBarry Smith 
1536*70c3da92SBarry Smith   /*  SuperLU is not currently supported through PETSc */
1537*70c3da92SBarry Smith #if defined(HAVE_SUPERLU)
1538*70c3da92SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
1539*70c3da92SBarry Smith   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
1540*70c3da92SBarry Smith #endif
1541*70c3da92SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
1542*70c3da92SBarry Smith   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
1543*70c3da92SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
1544*70c3da92SBarry Smith   if (flg) {
1545*70c3da92SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1546*70c3da92SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
1547*70c3da92SBarry Smith   }
1548*70c3da92SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1549*70c3da92SBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1550*70c3da92SBarry Smith   return 0;
1551*70c3da92SBarry Smith }
1552*70c3da92SBarry Smith 
1553*70c3da92SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
1554*70c3da92SBarry Smith {
1555*70c3da92SBarry Smith   Mat        C;
1556*70c3da92SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1557*70c3da92SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
1558*70c3da92SBarry Smith 
1559*70c3da92SBarry Smith   *B = 0;
1560*70c3da92SBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1561*70c3da92SBarry Smith   PLogObjectCreate(C);
1562*70c3da92SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1563*70c3da92SBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1564*70c3da92SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1565*70c3da92SBarry Smith   C->view       = MatView_SeqAIJ;
1566*70c3da92SBarry Smith   C->factor     = A->factor;
1567*70c3da92SBarry Smith   c->row        = 0;
1568*70c3da92SBarry Smith   c->col        = 0;
1569*70c3da92SBarry Smith   c->indexshift = shift;
1570*70c3da92SBarry Smith   C->assembled  = PETSC_TRUE;
1571*70c3da92SBarry Smith 
1572*70c3da92SBarry Smith   c->m = C->m   = a->m;
1573*70c3da92SBarry Smith   c->n = C->n   = a->n;
1574*70c3da92SBarry Smith   C->M          = a->m;
1575*70c3da92SBarry Smith   C->N          = a->n;
1576*70c3da92SBarry Smith 
1577*70c3da92SBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1578*70c3da92SBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
1579*70c3da92SBarry Smith   for ( i=0; i<m; i++ ) {
1580*70c3da92SBarry Smith     c->imax[i] = a->imax[i];
1581*70c3da92SBarry Smith     c->ilen[i] = a->ilen[i];
1582*70c3da92SBarry Smith   }
1583*70c3da92SBarry Smith 
1584*70c3da92SBarry Smith   /* allocate the matrix space */
1585*70c3da92SBarry Smith   c->singlemalloc = 1;
1586*70c3da92SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1587*70c3da92SBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1588*70c3da92SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1589*70c3da92SBarry Smith   c->i  = c->j + a->i[m] + shift;
1590*70c3da92SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
1591*70c3da92SBarry Smith   if (m > 0) {
1592*70c3da92SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1593*70c3da92SBarry Smith     if (cpvalues == COPY_VALUES) {
1594*70c3da92SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
1595*70c3da92SBarry Smith     }
1596*70c3da92SBarry Smith   }
1597*70c3da92SBarry Smith 
1598*70c3da92SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1599*70c3da92SBarry Smith   c->sorted      = a->sorted;
1600*70c3da92SBarry Smith   c->roworiented = a->roworiented;
1601*70c3da92SBarry Smith   c->nonew       = a->nonew;
1602*70c3da92SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
1603*70c3da92SBarry Smith 
1604*70c3da92SBarry Smith   if (a->diag) {
1605*70c3da92SBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1606*70c3da92SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
1607*70c3da92SBarry Smith     for ( i=0; i<m; i++ ) {
1608*70c3da92SBarry Smith       c->diag[i] = a->diag[i];
1609*70c3da92SBarry Smith     }
1610*70c3da92SBarry Smith   }
1611*70c3da92SBarry Smith   else c->diag          = 0;
1612*70c3da92SBarry Smith   c->inode.limit        = a->inode.limit;
1613*70c3da92SBarry Smith   c->inode.max_limit    = a->inode.max_limit;
1614*70c3da92SBarry Smith   if (a->inode.size){
1615*70c3da92SBarry Smith     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1616*70c3da92SBarry Smith     c->inode.node_count = a->inode.node_count;
1617*70c3da92SBarry Smith     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1618*70c3da92SBarry Smith   } else {
1619*70c3da92SBarry Smith     c->inode.size       = 0;
1620*70c3da92SBarry Smith     c->inode.node_count = 0;
1621*70c3da92SBarry Smith   }
1622*70c3da92SBarry Smith   c->nz                 = a->nz;
1623*70c3da92SBarry Smith   c->maxnz              = a->maxnz;
1624*70c3da92SBarry Smith   c->solve_work         = 0;
1625*70c3da92SBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1626*70c3da92SBarry Smith 
1627*70c3da92SBarry Smith   *B = C;
1628*70c3da92SBarry Smith   return 0;
1629*70c3da92SBarry Smith }
1630*70c3da92SBarry Smith 
1631*70c3da92SBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
1632*70c3da92SBarry Smith {
1633*70c3da92SBarry Smith   Mat_SeqAIJ   *a;
1634*70c3da92SBarry Smith   Mat          B;
1635*70c3da92SBarry Smith   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1636*70c3da92SBarry Smith   MPI_Comm     comm;
1637*70c3da92SBarry Smith 
1638*70c3da92SBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
1639*70c3da92SBarry Smith   MPI_Comm_size(comm,&size);
1640*70c3da92SBarry Smith   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
1641*70c3da92SBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
1642*70c3da92SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1643*70c3da92SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
1644*70c3da92SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
1645*70c3da92SBarry Smith 
1646*70c3da92SBarry Smith   /* read in row lengths */
1647*70c3da92SBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1648*70c3da92SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
1649*70c3da92SBarry Smith 
1650*70c3da92SBarry Smith   /* create our matrix */
1651*70c3da92SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1652*70c3da92SBarry Smith   B = *A;
1653*70c3da92SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1654*70c3da92SBarry Smith   shift = a->indexshift;
1655*70c3da92SBarry Smith 
1656*70c3da92SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1657*70c3da92SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
1658*70c3da92SBarry Smith   if (shift) {
1659*70c3da92SBarry Smith     for ( i=0; i<nz; i++ ) {
1660*70c3da92SBarry Smith       a->j[i] += 1;
1661*70c3da92SBarry Smith     }
1662*70c3da92SBarry Smith   }
1663*70c3da92SBarry Smith 
1664*70c3da92SBarry Smith   /* read in nonzero values */
1665*70c3da92SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
1666*70c3da92SBarry Smith 
1667*70c3da92SBarry Smith   /* set matrix "i" values */
1668*70c3da92SBarry Smith   a->i[0] = -shift;
1669*70c3da92SBarry Smith   for ( i=1; i<= M; i++ ) {
1670*70c3da92SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1671*70c3da92SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
1672*70c3da92SBarry Smith   }
1673*70c3da92SBarry Smith   PetscFree(rowlengths);
1674*70c3da92SBarry Smith 
1675*70c3da92SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1676*70c3da92SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1677*70c3da92SBarry Smith   return 0;
1678*70c3da92SBarry Smith }
1679*70c3da92SBarry Smith 
1680*70c3da92SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
1681*70c3da92SBarry Smith {
1682*70c3da92SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
1683*70c3da92SBarry Smith 
1684*70c3da92SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type");
1685*70c3da92SBarry Smith 
1686*70c3da92SBarry Smith   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
1687*70c3da92SBarry Smith   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1688*70c3da92SBarry Smith       (a->indexshift != b->indexshift)) {
1689*70c3da92SBarry Smith     *flg = PETSC_FALSE; return 0;
1690*70c3da92SBarry Smith   }
1691*70c3da92SBarry Smith 
1692*70c3da92SBarry Smith   /* if the a->i are the same */
1693*70c3da92SBarry Smith   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
1694*70c3da92SBarry Smith     *flg = PETSC_FALSE; return 0;
1695*70c3da92SBarry Smith   }
1696*70c3da92SBarry Smith 
1697*70c3da92SBarry Smith   /* if a->j are the same */
1698*70c3da92SBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
1699*70c3da92SBarry Smith     *flg = PETSC_FALSE; return 0;
1700*70c3da92SBarry Smith   }
1701*70c3da92SBarry Smith 
1702*70c3da92SBarry Smith   /* if a->a are the same */
1703*70c3da92SBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
1704*70c3da92SBarry Smith     *flg = PETSC_FALSE; return 0;
1705*70c3da92SBarry Smith   }
1706*70c3da92SBarry Smith   *flg = PETSC_TRUE;
1707*70c3da92SBarry Smith   return 0;
1708*70c3da92SBarry Smith 
1709*70c3da92SBarry Smith }
1710