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