xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 6d84be18fbb99ba69be7b8bdde5411a66955b7ea)
1*6d84be18SBarry Smith 
217ab2063SBarry Smith #ifndef lint
3*6d84be18SBarry Smith static char vcid[] = "$Id: aij.c,v 1.152 1996/03/01 03:43:59 bsmith Exp bsmith $";
417ab2063SBarry Smith #endif
517ab2063SBarry Smith 
6d5d45c9bSBarry Smith /*
7d5d45c9bSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
8d5d45c9bSBarry Smith   matrix storage format.
9d5d45c9bSBarry Smith */
1017ab2063SBarry Smith #include "aij.h"
1117ab2063SBarry Smith #include "vec/vecimpl.h"
1217ab2063SBarry Smith #include "inline/spops.h"
13e4d965acSSatish Balay #include "petsc.h"
1406763907SSatish Balay #include "inline/bitarray.h"
1517ab2063SBarry Smith 
167a743949SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int**,int**);
1717ab2063SBarry Smith 
18416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
1917ab2063SBarry Smith {
20416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
21a2744918SBarry Smith   int        ierr, *ia, *ja,n,*idx,i;
223d54f168SSatish Balay   /*Viewer     V1, V2;*/
2317ab2063SBarry Smith 
24a2744918SBarry Smith   /*
25a2744918SBarry Smith      this is tacky: In the future when we have written special factorization
26a2744918SBarry Smith      and solve routines for the identity permutation we should use a
27a2744918SBarry Smith      stride index set instead of the general one.
28a2744918SBarry Smith   */
29a2744918SBarry Smith   if (type  == ORDER_NATURAL) {
30a2744918SBarry Smith     n = a->n;
31a2744918SBarry Smith     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
32a2744918SBarry Smith     for ( i=0; i<n; i++ ) idx[i] = i;
33a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
34a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
35a2744918SBarry Smith     PetscFree(idx);
36a2744918SBarry Smith     ISSetPermutation(*rperm);
37a2744918SBarry Smith     ISSetPermutation(*cperm);
38a2744918SBarry Smith     ISSetIdentity(*rperm);
39a2744918SBarry Smith     ISSetIdentity(*cperm);
40a2744918SBarry Smith     return 0;
41a2744918SBarry Smith   }
42a2744918SBarry Smith 
437a743949SBarry Smith   ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,a->indexshift, &ia, &ja);CHKERRQ(ierr);
44416022c9SBarry Smith   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
45f31bba2fSSatish Balay 
460452661fSBarry Smith   PetscFree(ia); PetscFree(ja);
4717ab2063SBarry Smith   return 0;
4817ab2063SBarry Smith }
4917ab2063SBarry Smith 
50227d817aSBarry Smith #define CHUNKSIZE   15
5117ab2063SBarry Smith 
5217ab2063SBarry Smith /* This version has row oriented v  */
53416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
5417ab2063SBarry Smith {
55416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
56416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
574b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
58d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
59416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
6017ab2063SBarry Smith 
6117ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
62416022c9SBarry Smith     row  = im[k];
6317ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
64416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
6517ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
6617ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
67416022c9SBarry Smith     low = 0;
6817ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
69416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
70416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
714b0e389bSBarry Smith       col = in[l] - shift;
724b0e389bSBarry Smith       if (roworiented) {
734b0e389bSBarry Smith         value = *v++;
744b0e389bSBarry Smith       }
754b0e389bSBarry Smith       else {
764b0e389bSBarry Smith         value = v[k + l*m];
774b0e389bSBarry Smith       }
78416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
79416022c9SBarry Smith       while (high-low > 5) {
80416022c9SBarry Smith         t = (low+high)/2;
81416022c9SBarry Smith         if (rp[t] > col) high = t;
82416022c9SBarry Smith         else             low  = t;
8317ab2063SBarry Smith       }
84416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
8517ab2063SBarry Smith         if (rp[i] > col) break;
8617ab2063SBarry Smith         if (rp[i] == col) {
87416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
8817ab2063SBarry Smith           else                  ap[i] = value;
8917ab2063SBarry Smith           goto noinsert;
9017ab2063SBarry Smith         }
9117ab2063SBarry Smith       }
9217ab2063SBarry Smith       if (nonew) goto noinsert;
9317ab2063SBarry Smith       if (nrow >= rmax) {
9417ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
95416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
9617ab2063SBarry Smith         Scalar *new_a;
9717ab2063SBarry Smith 
9817ab2063SBarry Smith         /* malloc new storage space */
99416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1000452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
10117ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
10217ab2063SBarry Smith         new_i   = new_j + new_nz;
10317ab2063SBarry Smith 
10417ab2063SBarry Smith         /* copy over old data into new slots */
10517ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
106416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
107416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
108416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
109416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
11017ab2063SBarry Smith                                                            len*sizeof(int));
111416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
112416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
11317ab2063SBarry Smith                                                            len*sizeof(Scalar));
11417ab2063SBarry Smith         /* free up old matrix storage */
1150452661fSBarry Smith         PetscFree(a->a);
1160452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
117416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
118416022c9SBarry Smith         a->singlemalloc = 1;
11917ab2063SBarry Smith 
12017ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
121416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
122416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
123416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
124b810aeb4SBarry Smith         a->reallocs++;
12517ab2063SBarry Smith       }
126416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
127416022c9SBarry Smith       /* shift up all the later entries in this row */
128416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
12917ab2063SBarry Smith         rp[ii+1] = rp[ii];
13017ab2063SBarry Smith         ap[ii+1] = ap[ii];
13117ab2063SBarry Smith       }
13217ab2063SBarry Smith       rp[i] = col;
13317ab2063SBarry Smith       ap[i] = value;
13417ab2063SBarry Smith       noinsert:;
135416022c9SBarry Smith       low = i + 1;
13617ab2063SBarry Smith     }
13717ab2063SBarry Smith     ailen[row] = nrow;
13817ab2063SBarry Smith   }
13917ab2063SBarry Smith   return 0;
14017ab2063SBarry Smith }
14117ab2063SBarry Smith 
1427eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1437eb43aa7SLois Curfman McInnes {
1447eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
145b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1467eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1477eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
1487eb43aa7SLois Curfman McInnes 
1497eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
1507eb43aa7SLois Curfman McInnes     row  = im[k];
1517eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
1527eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
1537eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
1547eb43aa7SLois Curfman McInnes     nrow = ailen[row];
1557eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
1567eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
1577eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
1587eb43aa7SLois Curfman McInnes       col = in[l] - shift;
1597eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
1607eb43aa7SLois Curfman McInnes       while (high-low > 5) {
1617eb43aa7SLois Curfman McInnes         t = (low+high)/2;
1627eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
1637eb43aa7SLois Curfman McInnes         else             low  = t;
1647eb43aa7SLois Curfman McInnes       }
1657eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
1667eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
1677eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
168b49de8d1SLois Curfman McInnes           *v++ = ap[i];
1697eb43aa7SLois Curfman McInnes           goto finished;
1707eb43aa7SLois Curfman McInnes         }
1717eb43aa7SLois Curfman McInnes       }
172b49de8d1SLois Curfman McInnes       *v++ = zero;
1737eb43aa7SLois Curfman McInnes       finished:;
1747eb43aa7SLois Curfman McInnes     }
1757eb43aa7SLois Curfman McInnes   }
1767eb43aa7SLois Curfman McInnes   return 0;
1777eb43aa7SLois Curfman McInnes }
1787eb43aa7SLois Curfman McInnes 
17917ab2063SBarry Smith #include "draw.h"
18017ab2063SBarry Smith #include "pinclude/pviewer.h"
181416022c9SBarry Smith #include "sysio.h"
18217ab2063SBarry Smith 
183416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
18417ab2063SBarry Smith {
185416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
186416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
18717ab2063SBarry Smith 
188416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1890452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
190416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
191416022c9SBarry Smith   col_lens[1] = a->m;
192416022c9SBarry Smith   col_lens[2] = a->n;
193416022c9SBarry Smith   col_lens[3] = a->nz;
194416022c9SBarry Smith 
195416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
196416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
197416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
19817ab2063SBarry Smith   }
199416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
2000452661fSBarry Smith   PetscFree(col_lens);
201416022c9SBarry Smith 
202416022c9SBarry Smith   /* store column indices (zero start index) */
203416022c9SBarry Smith   if (a->indexshift) {
204416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
20517ab2063SBarry Smith   }
206416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
207416022c9SBarry Smith   if (a->indexshift) {
208416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
20917ab2063SBarry Smith   }
210416022c9SBarry Smith 
211416022c9SBarry Smith   /* store nonzero values */
212416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
21317ab2063SBarry Smith   return 0;
21417ab2063SBarry Smith }
215416022c9SBarry Smith 
216416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
217416022c9SBarry Smith {
218416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
21995e01e2fSLois Curfman McInnes   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg;
22017ab2063SBarry Smith   FILE        *fd;
22117ab2063SBarry Smith   char        *outputname;
22217ab2063SBarry Smith 
223227d817aSBarry Smith   ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr);
224416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
225416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
22617ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
22795e01e2fSLois Curfman McInnes     return 0;
22895e01e2fSLois Curfman McInnes   }
22995e01e2fSLois Curfman McInnes   else if (format == FILE_FORMAT_INFO_DETAILED) {
23095e01e2fSLois Curfman McInnes     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr);
23195e01e2fSLois Curfman McInnes     if (flg) fprintf(fd,"  not using I-node routines\n");
23295e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
23395e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
23417ab2063SBarry Smith   }
23517ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
23617ab2063SBarry Smith     int nz, nzalloc, mem;
237416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
238416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
23917ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
24017ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
24117ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
24217ab2063SBarry Smith 
24317ab2063SBarry Smith     for (i=0; i<m; i++) {
244416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
24517ab2063SBarry Smith #if defined(PETSC_COMPLEX)
2467a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]),
247416022c9SBarry Smith                    imag(a->a[j]));
24817ab2063SBarry Smith #else
2497a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
25017ab2063SBarry Smith #endif
25117ab2063SBarry Smith       }
25217ab2063SBarry Smith     }
25317ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
25417ab2063SBarry Smith   }
25517ab2063SBarry Smith   else {
25617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
25717ab2063SBarry Smith       fprintf(fd,"row %d:",i);
258416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
25917ab2063SBarry Smith #if defined(PETSC_COMPLEX)
260416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
261416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
26217ab2063SBarry Smith         }
26317ab2063SBarry Smith         else {
264416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
26517ab2063SBarry Smith         }
26617ab2063SBarry Smith #else
267416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
26817ab2063SBarry Smith #endif
26917ab2063SBarry Smith       }
27017ab2063SBarry Smith       fprintf(fd,"\n");
27117ab2063SBarry Smith     }
27217ab2063SBarry Smith   }
27317ab2063SBarry Smith   fflush(fd);
274416022c9SBarry Smith   return 0;
275416022c9SBarry Smith }
276416022c9SBarry Smith 
277416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
278416022c9SBarry Smith {
279416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
280cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
281cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
282d7e8b826SBarry Smith   Draw     draw = (Draw) viewer;
283cddf8d76SBarry Smith   DrawButton  button;
284cddf8d76SBarry Smith 
285416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
286416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
287416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
288416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
289cddf8d76SBarry Smith   color = DRAW_BLUE;
290416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
291cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
292416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
293cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
294cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
295cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
296cddf8d76SBarry Smith #else
297cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
298cddf8d76SBarry Smith #endif
299cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
300cddf8d76SBarry Smith     }
301cddf8d76SBarry Smith   }
302cddf8d76SBarry Smith   color = DRAW_CYAN;
303cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
304cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
305cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
306cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
307cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
308cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
309cddf8d76SBarry Smith     }
310cddf8d76SBarry Smith   }
311cddf8d76SBarry Smith   color = DRAW_RED;
312cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
313cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
314cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
315cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
316cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
317cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
318cddf8d76SBarry Smith #else
319cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
320cddf8d76SBarry Smith #endif
321cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
322416022c9SBarry Smith     }
323416022c9SBarry Smith   }
324416022c9SBarry Smith   DrawFlush(draw);
325cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
326cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
327cddf8d76SBarry Smith 
328cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
329cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
330cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
331cddf8d76SBarry Smith     DrawClear(draw);
332cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
333cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
334cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
335cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
336cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
337cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
338cddf8d76SBarry Smith     w *= scale; h *= scale;
339cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
340cddf8d76SBarry Smith     color = DRAW_BLUE;
341cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
342cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
343cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
344cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
345cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
346cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
347cddf8d76SBarry Smith #else
348cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
349cddf8d76SBarry Smith #endif
350cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
351cddf8d76SBarry Smith       }
352cddf8d76SBarry Smith     }
353cddf8d76SBarry Smith     color = DRAW_CYAN;
354cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
355cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
356cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
357cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
358cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
359cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
360cddf8d76SBarry Smith       }
361cddf8d76SBarry Smith     }
362cddf8d76SBarry Smith     color = DRAW_RED;
363cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
364cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
365cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
366cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
367cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
368cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
369cddf8d76SBarry Smith #else
370cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
371cddf8d76SBarry Smith #endif
372cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
373cddf8d76SBarry Smith       }
374cddf8d76SBarry Smith     }
375cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
376cddf8d76SBarry Smith   }
377416022c9SBarry Smith   return 0;
378416022c9SBarry Smith }
379416022c9SBarry Smith 
380416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
381416022c9SBarry Smith {
382416022c9SBarry Smith   Mat         A = (Mat) obj;
383416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
384416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
385416022c9SBarry Smith 
386416022c9SBarry Smith   if (!viewer) {
387416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
388416022c9SBarry Smith   }
389416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
390416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
391416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
392416022c9SBarry Smith     }
393416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
394416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
395416022c9SBarry Smith     }
396416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
397416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
398416022c9SBarry Smith     }
399416022c9SBarry Smith   }
400416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
401416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
402416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
40317ab2063SBarry Smith   }
40417ab2063SBarry Smith   return 0;
40517ab2063SBarry Smith }
406c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
407416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
40817ab2063SBarry Smith {
409416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
41041c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
411416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
412416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
41317ab2063SBarry Smith 
41417ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
41517ab2063SBarry Smith 
41617ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
417416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
41817ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
41917ab2063SBarry Smith     if (fshift) {
420416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
42117ab2063SBarry Smith       N = ailen[i];
42217ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
42317ab2063SBarry Smith         ip[j-fshift] = ip[j];
42417ab2063SBarry Smith         ap[j-fshift] = ap[j];
42517ab2063SBarry Smith       }
42617ab2063SBarry Smith     }
42717ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
42817ab2063SBarry Smith   }
42917ab2063SBarry Smith   if (m) {
43017ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
43117ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
43217ab2063SBarry Smith   }
43317ab2063SBarry Smith   /* reset ilen and imax for each row */
43417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
43517ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
43617ab2063SBarry Smith   }
437416022c9SBarry Smith   a->nz = ai[m] + shift;
43817ab2063SBarry Smith 
43917ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
440416022c9SBarry Smith   if (fshift && a->diag) {
4410452661fSBarry Smith     PetscFree(a->diag);
442416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
443416022c9SBarry Smith     a->diag = 0;
44417ab2063SBarry Smith   }
445b810aeb4SBarry Smith   PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n",
446b810aeb4SBarry Smith            fshift,a->nz,m);
447b810aeb4SBarry Smith   PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n",
448b810aeb4SBarry Smith            a->reallocs);
44976dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
45041c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
45117ab2063SBarry Smith   return 0;
45217ab2063SBarry Smith }
45317ab2063SBarry Smith 
454416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
45517ab2063SBarry Smith {
456416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
457cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
45817ab2063SBarry Smith   return 0;
45917ab2063SBarry Smith }
460416022c9SBarry Smith 
46117ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
46217ab2063SBarry Smith {
463416022c9SBarry Smith   Mat        A  = (Mat) obj;
464416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
465d5d45c9bSBarry Smith 
46617ab2063SBarry Smith #if defined(PETSC_LOG)
467416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
46817ab2063SBarry Smith #endif
4690452661fSBarry Smith   PetscFree(a->a);
4700452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4710452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
4720452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4730452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
4740452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
47576dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
4760452661fSBarry Smith   PetscFree(a);
477416022c9SBarry Smith   PLogObjectDestroy(A);
4780452661fSBarry Smith   PetscHeaderDestroy(A);
47917ab2063SBarry Smith   return 0;
48017ab2063SBarry Smith }
48117ab2063SBarry Smith 
482416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
48317ab2063SBarry Smith {
48417ab2063SBarry Smith   return 0;
48517ab2063SBarry Smith }
48617ab2063SBarry Smith 
487416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
48817ab2063SBarry Smith {
489416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
490416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
4914b0e389bSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
492416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
493416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
494416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
495e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
496e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
497e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
498e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
499df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
500df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
5014d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
5026c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
5036c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
5046c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
5056c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
5066c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
507e2f28af5SBarry Smith   else
5084d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
50917ab2063SBarry Smith   return 0;
51017ab2063SBarry Smith }
51117ab2063SBarry Smith 
512416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
51317ab2063SBarry Smith {
514416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
515416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
51617ab2063SBarry Smith   Scalar     *x, zero = 0.0;
51717ab2063SBarry Smith 
51817ab2063SBarry Smith   VecSet(&zero,v);
51917ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
520416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
521416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
522416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
523416022c9SBarry Smith       if (a->j[j]+shift == i) {
524416022c9SBarry Smith         x[i] = a->a[j];
52517ab2063SBarry Smith         break;
52617ab2063SBarry Smith       }
52717ab2063SBarry Smith     }
52817ab2063SBarry Smith   }
52917ab2063SBarry Smith   return 0;
53017ab2063SBarry Smith }
53117ab2063SBarry Smith 
53217ab2063SBarry Smith /* -------------------------------------------------------*/
53317ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
53417ab2063SBarry Smith /* -------------------------------------------------------*/
535416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
53617ab2063SBarry Smith {
537416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
53817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
539416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
54017ab2063SBarry Smith 
54117ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
542cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
54317ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
54417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
545416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
546416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
547416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
54817ab2063SBarry Smith     alpha = x[i];
54917ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
55017ab2063SBarry Smith   }
551416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
55217ab2063SBarry Smith   return 0;
55317ab2063SBarry Smith }
554d5d45c9bSBarry Smith 
555416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
55617ab2063SBarry Smith {
557416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
55817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
559416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
56017ab2063SBarry Smith 
56117ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
56217ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
56317ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
56417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
565416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
566416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
567416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
56817ab2063SBarry Smith     alpha = x[i];
56917ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
57017ab2063SBarry Smith   }
57117ab2063SBarry Smith   return 0;
57217ab2063SBarry Smith }
57317ab2063SBarry Smith 
574416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
57517ab2063SBarry Smith {
576416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
57717ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
578416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
57917ab2063SBarry Smith 
58017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
58117ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
582416022c9SBarry Smith   idx  = a->j;
583416022c9SBarry Smith   v    = a->a;
584416022c9SBarry Smith   ii   = a->i;
58517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
586416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
58717ab2063SBarry Smith     sum  = 0.0;
58817ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
58917ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
590416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
59117ab2063SBarry Smith     y[i] = sum;
59217ab2063SBarry Smith   }
593416022c9SBarry Smith   PLogFlops(2*a->nz - m);
59417ab2063SBarry Smith   return 0;
59517ab2063SBarry Smith }
59617ab2063SBarry Smith 
597416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
59817ab2063SBarry Smith {
599416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
60017ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
601cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
60217ab2063SBarry Smith 
60317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
60417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
605cddf8d76SBarry Smith   idx  = a->j;
606cddf8d76SBarry Smith   v    = a->a;
607cddf8d76SBarry Smith   ii   = a->i;
60817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
609cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
61017ab2063SBarry Smith     sum  = y[i];
611cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
612cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
61317ab2063SBarry Smith     z[i] = sum;
61417ab2063SBarry Smith   }
615416022c9SBarry Smith   PLogFlops(2*a->nz);
61617ab2063SBarry Smith   return 0;
61717ab2063SBarry Smith }
61817ab2063SBarry Smith 
61917ab2063SBarry Smith /*
62017ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
62117ab2063SBarry Smith */
62217ab2063SBarry Smith 
623416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
62417ab2063SBarry Smith {
625416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
626416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
62717ab2063SBarry Smith 
6280452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
629416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
630416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
631416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
632416022c9SBarry Smith       if (a->j[j]+shift == i) {
63317ab2063SBarry Smith         diag[i] = j - shift;
63417ab2063SBarry Smith         break;
63517ab2063SBarry Smith       }
63617ab2063SBarry Smith     }
63717ab2063SBarry Smith   }
638416022c9SBarry Smith   a->diag = diag;
63917ab2063SBarry Smith   return 0;
64017ab2063SBarry Smith }
64117ab2063SBarry Smith 
642416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
64317ab2063SBarry Smith                            double fshift,int its,Vec xx)
64417ab2063SBarry Smith {
645416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
646416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
647d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
64817ab2063SBarry Smith 
64917ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
650416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
651416022c9SBarry Smith   diag = a->diag;
652416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
65317ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
65417ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
65517ab2063SBarry Smith     bs = b + shift;
65617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
657416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
658416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
659416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
660416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
66117ab2063SBarry Smith         sum  = b[i]*d/omega;
66217ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
66317ab2063SBarry Smith         x[i] = sum;
66417ab2063SBarry Smith     }
66517ab2063SBarry Smith     return 0;
66617ab2063SBarry Smith   }
66717ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
66817ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
66917ab2063SBarry Smith   }
670416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
67117ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
67217ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
67317ab2063SBarry Smith 
67417ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
67517ab2063SBarry Smith 
67617ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
67717ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
67817ab2063SBarry Smith     is the relaxation factor.
67917ab2063SBarry Smith     */
6800452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
68117ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
68217ab2063SBarry Smith 
68317ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
68417ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
685416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
686416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
687416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
688416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
68917ab2063SBarry Smith       sum  = b[i];
69017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
69117ab2063SBarry Smith       x[i] = omega*(sum/d);
69217ab2063SBarry Smith     }
69317ab2063SBarry Smith 
69417ab2063SBarry Smith     /*  t = b - (2*E - D)x */
695416022c9SBarry Smith     v = a->a;
69617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
69717ab2063SBarry Smith 
69817ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
699416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
700416022c9SBarry Smith     diag = a->diag;
70117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
702416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
703416022c9SBarry Smith       n    = diag[i] - a->i[i];
704416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
705416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
70617ab2063SBarry Smith       sum  = t[i];
70717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
70817ab2063SBarry Smith       t[i] = omega*(sum/d);
70917ab2063SBarry Smith     }
71017ab2063SBarry Smith 
71117ab2063SBarry Smith     /*  x = x + t */
71217ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7130452661fSBarry Smith     PetscFree(t);
71417ab2063SBarry Smith     return 0;
71517ab2063SBarry Smith   }
71617ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
71717ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
71817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
719416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
720416022c9SBarry Smith         n    = diag[i] - a->i[i];
721416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
722416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
72317ab2063SBarry Smith         sum  = b[i];
72417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
72517ab2063SBarry Smith         x[i] = omega*(sum/d);
72617ab2063SBarry Smith       }
72717ab2063SBarry Smith       xb = x;
72817ab2063SBarry Smith     }
72917ab2063SBarry Smith     else xb = b;
73017ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
73117ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
73217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
733416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
73417ab2063SBarry Smith       }
73517ab2063SBarry Smith     }
73617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
73717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
738416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
739416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
740416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
741416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
74217ab2063SBarry Smith         sum  = xb[i];
74317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
74417ab2063SBarry Smith         x[i] = omega*(sum/d);
74517ab2063SBarry Smith       }
74617ab2063SBarry Smith     }
74717ab2063SBarry Smith     its--;
74817ab2063SBarry Smith   }
74917ab2063SBarry Smith   while (its--) {
75017ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
75117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
752416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
753416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
754416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
755416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
75617ab2063SBarry Smith         sum  = b[i];
75717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
75817ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
75917ab2063SBarry Smith       }
76017ab2063SBarry Smith     }
76117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
76217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
763416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
764416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
765416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
766416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
76717ab2063SBarry Smith         sum  = b[i];
76817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
76917ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
77017ab2063SBarry Smith       }
77117ab2063SBarry Smith     }
77217ab2063SBarry Smith   }
77317ab2063SBarry Smith   return 0;
77417ab2063SBarry Smith }
77517ab2063SBarry Smith 
776d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
77717ab2063SBarry Smith {
778416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
779416022c9SBarry Smith   *nz      = a->nz;
780416022c9SBarry Smith   *nzalloc = a->maxnz;
781416022c9SBarry Smith   *mem     = (int)A->mem;
78217ab2063SBarry Smith   return 0;
78317ab2063SBarry Smith }
78417ab2063SBarry Smith 
78517ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
78617ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
78717ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
78817ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
78917ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
79017ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
79117ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
79217ab2063SBarry Smith 
79317ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
79417ab2063SBarry Smith {
795416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
796416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
79717ab2063SBarry Smith 
79817ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
79917ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
80017ab2063SBarry Smith   if (diag) {
80117ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
802416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
803416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
804416022c9SBarry Smith         a->ilen[rows[i]] = 1;
805416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
806416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
80717ab2063SBarry Smith       }
80817ab2063SBarry Smith       else {
80917ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
81017ab2063SBarry Smith         CHKERRQ(ierr);
81117ab2063SBarry Smith       }
81217ab2063SBarry Smith     }
81317ab2063SBarry Smith   }
81417ab2063SBarry Smith   else {
81517ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
816416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
817416022c9SBarry Smith       a->ilen[rows[i]] = 0;
81817ab2063SBarry Smith     }
81917ab2063SBarry Smith   }
82017ab2063SBarry Smith   ISRestoreIndices(is,&rows);
82117ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82217ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82317ab2063SBarry Smith   return 0;
82417ab2063SBarry Smith }
82517ab2063SBarry Smith 
826416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
82717ab2063SBarry Smith {
828416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
829416022c9SBarry Smith   *m = a->m; *n = a->n;
83017ab2063SBarry Smith   return 0;
83117ab2063SBarry Smith }
83217ab2063SBarry Smith 
833416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
83417ab2063SBarry Smith {
835416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
836416022c9SBarry Smith   *m = 0; *n = a->m;
83717ab2063SBarry Smith   return 0;
83817ab2063SBarry Smith }
8394e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
84017ab2063SBarry Smith {
841416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
842c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
84317ab2063SBarry Smith 
844416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
84517ab2063SBarry Smith 
846416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
847416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
84817ab2063SBarry Smith   if (idx) {
849416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
8504e093b46SBarry Smith     if (*nz && shift) {
8510452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
85217ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
8534e093b46SBarry Smith     } else if (*nz) {
8544e093b46SBarry Smith       *idx = itmp;
85517ab2063SBarry Smith     }
85617ab2063SBarry Smith     else *idx = 0;
85717ab2063SBarry Smith   }
85817ab2063SBarry Smith   return 0;
85917ab2063SBarry Smith }
86017ab2063SBarry Smith 
8614e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
86217ab2063SBarry Smith {
8634e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
8644e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
86517ab2063SBarry Smith   return 0;
86617ab2063SBarry Smith }
86717ab2063SBarry Smith 
868cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
86917ab2063SBarry Smith {
870416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
871416022c9SBarry Smith   Scalar     *v = a->a;
87217ab2063SBarry Smith   double     sum = 0.0;
873416022c9SBarry Smith   int        i, j,shift = a->indexshift;
87417ab2063SBarry Smith 
87517ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
876416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
87717ab2063SBarry Smith #if defined(PETSC_COMPLEX)
87817ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
87917ab2063SBarry Smith #else
88017ab2063SBarry Smith       sum += (*v)*(*v); v++;
88117ab2063SBarry Smith #endif
88217ab2063SBarry Smith     }
88317ab2063SBarry Smith     *norm = sqrt(sum);
88417ab2063SBarry Smith   }
88517ab2063SBarry Smith   else if (type == NORM_1) {
88617ab2063SBarry Smith     double *tmp;
887416022c9SBarry Smith     int    *jj = a->j;
8880452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
889cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
89017ab2063SBarry Smith     *norm = 0.0;
891416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
892a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
89317ab2063SBarry Smith     }
894416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
89517ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
89617ab2063SBarry Smith     }
8970452661fSBarry Smith     PetscFree(tmp);
89817ab2063SBarry Smith   }
89917ab2063SBarry Smith   else if (type == NORM_INFINITY) {
90017ab2063SBarry Smith     *norm = 0.0;
901416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
902416022c9SBarry Smith       v = a->a + a->i[j] + shift;
90317ab2063SBarry Smith       sum = 0.0;
904416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
905cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
90617ab2063SBarry Smith       }
90717ab2063SBarry Smith       if (sum > *norm) *norm = sum;
90817ab2063SBarry Smith     }
90917ab2063SBarry Smith   }
91017ab2063SBarry Smith   else {
91148d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
91217ab2063SBarry Smith   }
91317ab2063SBarry Smith   return 0;
91417ab2063SBarry Smith }
91517ab2063SBarry Smith 
916416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
91717ab2063SBarry Smith {
918416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
919416022c9SBarry Smith   Mat        C;
920416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
921416022c9SBarry Smith   int        shift = a->indexshift;
922d5d45c9bSBarry Smith   Scalar     *array = a->a;
92317ab2063SBarry Smith 
9243638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
9253638b69dSLois Curfman McInnes     SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place");
9260452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
927cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
92817ab2063SBarry Smith   if (shift) {
92917ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
93017ab2063SBarry Smith   }
93117ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
932416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
9330452661fSBarry Smith   PetscFree(col);
93417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
93517ab2063SBarry Smith     len = ai[i+1]-ai[i];
936416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
93717ab2063SBarry Smith     array += len; aj += len;
93817ab2063SBarry Smith   }
93917ab2063SBarry Smith   if (shift) {
94017ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
94117ab2063SBarry Smith   }
94217ab2063SBarry Smith 
943416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
944416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
94517ab2063SBarry Smith 
9463638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
947416022c9SBarry Smith     *B = C;
94817ab2063SBarry Smith   } else {
949416022c9SBarry Smith     /* This isn't really an in-place transpose */
9500452661fSBarry Smith     PetscFree(a->a);
9510452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
9520452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
9530452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
9540452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
9550452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
9560452661fSBarry Smith     PetscFree(a);
957416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
9580452661fSBarry Smith     PetscHeaderDestroy(C);
95917ab2063SBarry Smith   }
96017ab2063SBarry Smith   return 0;
96117ab2063SBarry Smith }
96217ab2063SBarry Smith 
963f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
96417ab2063SBarry Smith {
965416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
96617ab2063SBarry Smith   Scalar     *l,*r,x,*v;
967d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
96817ab2063SBarry Smith 
96917ab2063SBarry Smith   if (ll) {
97017ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
971f0b747eeSBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
972416022c9SBarry Smith     v = a->a;
97317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
97417ab2063SBarry Smith       x = l[i];
975416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
97617ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
97717ab2063SBarry Smith     }
97817ab2063SBarry Smith   }
97917ab2063SBarry Smith   if (rr) {
98017ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
981f0b747eeSBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
982416022c9SBarry Smith     v = a->a; jj = a->j;
98317ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
98417ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
98517ab2063SBarry Smith     }
98617ab2063SBarry Smith   }
98717ab2063SBarry Smith   return 0;
98817ab2063SBarry Smith }
98917ab2063SBarry Smith 
990cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
99117ab2063SBarry Smith {
992db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
99302834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
994a2744918SBarry Smith   register int sum,lensi;
99502834360SBarry Smith   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
996a2744918SBarry Smith   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
997db02288aSLois Curfman McInnes   Scalar       *vwork,*a_new;
998416022c9SBarry Smith   Mat          C;
99917ab2063SBarry Smith 
100017ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
100117ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
100217ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
100317ab2063SBarry Smith 
10047264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
100502834360SBarry Smith     /* special case of contiguous rows */
10060452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
100702834360SBarry Smith     starts = lens + ncols;
100802834360SBarry Smith     /* loop over new rows determining lens and starting points */
100902834360SBarry Smith     for (i=0; i<nrows; i++) {
1010a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1011a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
101202834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1013d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
101402834360SBarry Smith           starts[i] = k;
101502834360SBarry Smith           break;
101602834360SBarry Smith 	}
101702834360SBarry Smith       }
1018a2744918SBarry Smith       sum = 0;
101902834360SBarry Smith       while (k < kend) {
1020d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1021a2744918SBarry Smith         sum++;
102202834360SBarry Smith       }
1023a2744918SBarry Smith       lens[i] = sum;
102402834360SBarry Smith     }
102502834360SBarry Smith     /* create submatrix */
1026cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
102708480c60SBarry Smith       int n_cols,n_rows;
102808480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
102908480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1030d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
103108480c60SBarry Smith       C = *B;
103208480c60SBarry Smith     }
103308480c60SBarry Smith     else {
103402834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
103508480c60SBarry Smith     }
1036db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1037db02288aSLois Curfman McInnes 
103802834360SBarry Smith     /* loop over rows inserting into submatrix */
1039db02288aSLois Curfman McInnes     a_new    = c->a;
1040db02288aSLois Curfman McInnes     j_new    = c->j;
1041db02288aSLois Curfman McInnes     i_new    = c->i;
1042db02288aSLois Curfman McInnes     i_new[0] = -shift;
104302834360SBarry Smith     for (i=0; i<nrows; i++) {
1044a2744918SBarry Smith       ii    = starts[i];
1045a2744918SBarry Smith       lensi = lens[i];
1046a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1047a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
104802834360SBarry Smith       }
1049a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1050a2744918SBarry Smith       a_new      += lensi;
1051a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1052a2744918SBarry Smith       c->ilen[i]  = lensi;
105302834360SBarry Smith     }
10540452661fSBarry Smith     PetscFree(lens);
105502834360SBarry Smith   }
105602834360SBarry Smith   else {
105702834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10587264ac53SSatish Balay     ierr = SYIsort(ncols,icol); CHKERRQ(ierr);
10590452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
106002834360SBarry Smith     ssmap = smap + shift;
10617eb43aa7SLois Curfman McInnes     cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork);
106202834360SBarry Smith     lens  = cwork + ncols;
10630452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1064cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
106517ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
106602834360SBarry Smith     /* determine lens of each row */
106702834360SBarry Smith     for (i=0; i<nrows; i++) {
1068d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
106902834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
107002834360SBarry Smith       lens[i] = 0;
107102834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1072d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
107302834360SBarry Smith           lens[i]++;
107402834360SBarry Smith         }
107502834360SBarry Smith       }
107602834360SBarry Smith     }
107717ab2063SBarry Smith     /* Create and fill new matrix */
1078a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
107908480c60SBarry Smith       int n_cols,n_rows;
108008480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
108108480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1082d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
108308480c60SBarry Smith       C = *B;
108408480c60SBarry Smith     }
108508480c60SBarry Smith     else {
108602834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
108708480c60SBarry Smith     }
108817ab2063SBarry Smith     for (i=0; i<nrows; i++) {
108917ab2063SBarry Smith       nznew  = 0;
1090d8ced48eSBarry Smith       kstart = ai[irow[i]]+shift;
1091416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
109217ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
109302834360SBarry Smith         if (ssmap[a->j[k]]) {
109402834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1095416022c9SBarry Smith           vwork[nznew++] = a->a[k];
109617ab2063SBarry Smith         }
109717ab2063SBarry Smith       }
1098416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
109917ab2063SBarry Smith     }
110002834360SBarry Smith     /* Free work space */
110102834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
11020452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
110302834360SBarry Smith   }
1104416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1105416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
110617ab2063SBarry Smith 
110717ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1108416022c9SBarry Smith   *B = C;
110917ab2063SBarry Smith   return 0;
111017ab2063SBarry Smith }
111117ab2063SBarry Smith 
1112a871dcd8SBarry Smith /*
111363b91edcSBarry Smith      note: This can only work for identity for row and col. It would
111463b91edcSBarry Smith    be good to check this and otherwise generate an error.
1115a871dcd8SBarry Smith */
111663b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1117a871dcd8SBarry Smith {
111863b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
111908480c60SBarry Smith   int        ierr;
112063b91edcSBarry Smith   Mat        outA;
112163b91edcSBarry Smith 
1122a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1123a871dcd8SBarry Smith 
112463b91edcSBarry Smith   outA          = inA;
112563b91edcSBarry Smith   inA->factor   = FACTOR_LU;
112663b91edcSBarry Smith   a->row        = row;
112763b91edcSBarry Smith   a->col        = col;
112863b91edcSBarry Smith 
11290452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
113063b91edcSBarry Smith 
113108480c60SBarry Smith   if (!a->diag) {
113208480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
113363b91edcSBarry Smith   }
113463b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1135a871dcd8SBarry Smith   return 0;
1136a871dcd8SBarry Smith }
1137a871dcd8SBarry Smith 
1138f0b747eeSBarry Smith #include "pinclude/plapack.h"
1139f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1140f0b747eeSBarry Smith {
1141f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1142f0b747eeSBarry Smith   int        one = 1;
1143f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1144f0b747eeSBarry Smith   PLogFlops(a->nz);
1145f0b747eeSBarry Smith   return 0;
1146f0b747eeSBarry Smith }
1147f0b747eeSBarry Smith 
1148cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1149cddf8d76SBarry Smith                                     Mat **B)
1150cddf8d76SBarry Smith {
1151cddf8d76SBarry Smith   int ierr,i;
1152cddf8d76SBarry Smith 
1153cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
11540452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1155cddf8d76SBarry Smith   }
1156cddf8d76SBarry Smith 
1157cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1158cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1159cddf8d76SBarry Smith   }
1160cddf8d76SBarry Smith   return 0;
1161cddf8d76SBarry Smith }
1162cddf8d76SBarry Smith 
1163e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
11644dcbc457SBarry Smith {
1165e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
116606763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
11678a047759SSatish Balay   int        start, end, *ai, *aj;
116806763907SSatish Balay   char       *table;
11698a047759SSatish Balay   shift = a->indexshift;
1170e4d965acSSatish Balay   m     = a->m;
1171e4d965acSSatish Balay   ai    = a->i;
11728a047759SSatish Balay   aj    = a->j+shift;
11738a047759SSatish Balay 
1174e4d965acSSatish Balay   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
117506763907SSatish Balay 
117606763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
117706763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
117806763907SSatish Balay 
1179e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1180e4d965acSSatish Balay     /* Initialise the two local arrays */
1181e4d965acSSatish Balay     isz  = 0;
118206763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1183e4d965acSSatish Balay 
1184e4d965acSSatish Balay                 /* Extract the indices, assume there can be duplicate entries */
11854dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
1186e4d965acSSatish Balay     ierr = ISGetLocalSize(is[i],&n);  CHKERRQ(ierr);
1187e4d965acSSatish Balay 
1188e4d965acSSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1189e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
119006763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
11914dcbc457SBarry Smith     }
119206763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
119306763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1194e4d965acSSatish Balay 
119504a348a9SBarry Smith     k = 0;
119604a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
119704a348a9SBarry Smith       n = isz;
119806763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1199e4d965acSSatish Balay         row   = nidx[k];
1200e4d965acSSatish Balay         start = ai[row];
1201e4d965acSSatish Balay         end   = ai[row+1];
120204a348a9SBarry Smith         for ( l = start; l<end ; l++){
12038a047759SSatish Balay           val = aj[l] + shift;
120406763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1205e4d965acSSatish Balay         }
1206e4d965acSSatish Balay       }
1207e4d965acSSatish Balay     }
1208e4d965acSSatish Balay     ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1209e4d965acSSatish Balay   }
121004a348a9SBarry Smith   PetscFree(table);
121106763907SSatish Balay   PetscFree(nidx);
1212e4d965acSSatish Balay   return 0;
12134dcbc457SBarry Smith }
121417ab2063SBarry Smith 
1215682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1216682d7d0cSBarry Smith {
1217682d7d0cSBarry Smith   static int called = 0;
1218682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1219682d7d0cSBarry Smith 
1220682d7d0cSBarry Smith   if (called) return 0; else called = 1;
12212a7368beSLois Curfman McInnes   MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1222682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
12232a7368beSLois Curfman McInnes   MPIU_printf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
1224682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
1225682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1226682d7d0cSBarry Smith #if defined(HAVE_ESSL)
1227682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1228682d7d0cSBarry Smith #endif
1229682d7d0cSBarry Smith   return 0;
1230682d7d0cSBarry Smith }
12317264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg);
1232682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
123317ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
123417ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1235416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1236416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
123717ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
123817ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
123917ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
124017ab2063SBarry Smith        MatRelax_SeqAIJ,
124117ab2063SBarry Smith        MatTranspose_SeqAIJ,
12427264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1243f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
124417ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
124517ab2063SBarry Smith        MatCompress_SeqAIJ,
124617ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
124717ab2063SBarry Smith        MatGetReordering_SeqAIJ,
124817ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
124917ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
125017ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
125117ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1252416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
12533d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1254cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
12557eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1256682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1257f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
1258f0b747eeSBarry Smith        MatScale_SeqAIJ};
125917ab2063SBarry Smith 
126017ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
126117ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
126217ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
126317ab2063SBarry Smith 
126417ab2063SBarry Smith /*@C
1265682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
12660d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
12676e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
12680d15e28bSLois Curfman McInnes    (or nzz).  By setting these parameters accurately, performance can be
12690d15e28bSLois Curfman McInnes    increased by more than a factor of 50.
127017ab2063SBarry Smith 
127117ab2063SBarry Smith    Input Parameters:
127217ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
127317ab2063SBarry Smith .  m - number of rows
127417ab2063SBarry Smith .  n - number of columns
127517ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
127617ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
127717ab2063SBarry Smith 
127817ab2063SBarry Smith    Output Parameter:
1279416022c9SBarry Smith .  A - the matrix
128017ab2063SBarry Smith 
128117ab2063SBarry Smith    Notes:
128217ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
128317ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
12840002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
12850002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
128617ab2063SBarry Smith 
128717ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1288a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
12890d15e28bSLois Curfman McInnes    allocation.  For additional details, see the users manual chapter on
12900d15e28bSLois Curfman McInnes    matrices and the file $(PETSC_DIR)/Performance.
129117ab2063SBarry Smith 
1292682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1293682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1294682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
12956c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
12966c7ebb05SLois Curfman McInnes 
12976c7ebb05SLois Curfman McInnes    Options Database Keys:
12986c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
12996e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
13006e62573dSLois Curfman McInnes $        (max limit=5)
13016e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
13026e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
13036e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
130417ab2063SBarry Smith 
130517ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
130617ab2063SBarry Smith @*/
1307416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
130817ab2063SBarry Smith {
1309416022c9SBarry Smith   Mat        B;
1310416022c9SBarry Smith   Mat_SeqAIJ *b;
131169957df2SSatish Balay   int        i,len,ierr, flg;
1312d5d45c9bSBarry Smith 
1313416022c9SBarry Smith   *A      = 0;
13140452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1315416022c9SBarry Smith   PLogObjectCreate(B);
13160452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1317416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1318416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1319416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1320416022c9SBarry Smith   B->factor           = 0;
1321416022c9SBarry Smith   B->lupivotthreshold = 1.0;
13227a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
132369957df2SSatish Balay                           &flg); CHKERRQ(ierr);
13247a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
13257a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
13267a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1327416022c9SBarry Smith   b->row              = 0;
1328416022c9SBarry Smith   b->col              = 0;
1329416022c9SBarry Smith   b->indexshift       = 0;
1330b810aeb4SBarry Smith   b->reallocs         = 0;
133169957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
133269957df2SSatish Balay   if (flg) b->indexshift = -1;
133317ab2063SBarry Smith 
1334416022c9SBarry Smith   b->m       = m;
1335416022c9SBarry Smith   b->n       = n;
13360452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1337b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
13387b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
13397b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1340416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
134117ab2063SBarry Smith     nz = nz*m;
134217ab2063SBarry Smith   }
134317ab2063SBarry Smith   else {
134417ab2063SBarry Smith     nz = 0;
1345416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
134617ab2063SBarry Smith   }
134717ab2063SBarry Smith 
134817ab2063SBarry Smith   /* allocate the matrix space */
1349416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
13500452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1351416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1352cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1353416022c9SBarry Smith   b->i  = b->j + nz;
1354416022c9SBarry Smith   b->singlemalloc = 1;
135517ab2063SBarry Smith 
1356416022c9SBarry Smith   b->i[0] = -b->indexshift;
135717ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1358416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
135917ab2063SBarry Smith   }
136017ab2063SBarry Smith 
1361416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
13620452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1363416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1364416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
136517ab2063SBarry Smith 
1366416022c9SBarry Smith   b->nz               = 0;
1367416022c9SBarry Smith   b->maxnz            = nz;
1368416022c9SBarry Smith   b->sorted           = 0;
1369416022c9SBarry Smith   b->roworiented      = 1;
1370416022c9SBarry Smith   b->nonew            = 0;
1371416022c9SBarry Smith   b->diag             = 0;
1372416022c9SBarry Smith   b->solve_work       = 0;
137371bd300dSLois Curfman McInnes   b->spptr            = 0;
1374754ec7b1SSatish Balay   b->inode.node_count = 0;
1375754ec7b1SSatish Balay   b->inode.size       = 0;
13766c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
13776c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
137817ab2063SBarry Smith 
1379416022c9SBarry Smith   *A = B;
13804b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
13814b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
138269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
138369957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
13844b14c69eSBarry Smith #endif
138569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
138669957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
138769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
138869957df2SSatish Balay   if (flg) {
1389416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1390416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
139117ab2063SBarry Smith   }
139269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
139369957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
139417ab2063SBarry Smith   return 0;
139517ab2063SBarry Smith }
139617ab2063SBarry Smith 
13973d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
139817ab2063SBarry Smith {
1399416022c9SBarry Smith   Mat        C;
1400416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
140108480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
140217ab2063SBarry Smith 
14034043dd9cSLois Curfman McInnes   *B = 0;
14040452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1405416022c9SBarry Smith   PLogObjectCreate(C);
14060452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
140741c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1408416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1409416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1410416022c9SBarry Smith   C->factor     = A->factor;
1411416022c9SBarry Smith   c->row        = 0;
1412416022c9SBarry Smith   c->col        = 0;
1413416022c9SBarry Smith   c->indexshift = shift;
1414c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
141517ab2063SBarry Smith 
1416416022c9SBarry Smith   c->m          = a->m;
1417416022c9SBarry Smith   c->n          = a->n;
141817ab2063SBarry Smith 
14190452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
14200452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
142117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1422416022c9SBarry Smith     c->imax[i] = a->imax[i];
1423416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
142417ab2063SBarry Smith   }
142517ab2063SBarry Smith 
142617ab2063SBarry Smith   /* allocate the matrix space */
1427416022c9SBarry Smith   c->singlemalloc = 1;
1428416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
14290452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1430416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1431416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1432416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
143317ab2063SBarry Smith   if (m > 0) {
1434416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
143508480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1436416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
143717ab2063SBarry Smith     }
143808480c60SBarry Smith   }
143917ab2063SBarry Smith 
1440416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1441416022c9SBarry Smith   c->sorted      = a->sorted;
1442416022c9SBarry Smith   c->roworiented = a->roworiented;
1443416022c9SBarry Smith   c->nonew       = a->nonew;
14447a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
144517ab2063SBarry Smith 
1446416022c9SBarry Smith   if (a->diag) {
14470452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1448416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
144917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1450416022c9SBarry Smith       c->diag[i] = a->diag[i];
145117ab2063SBarry Smith     }
145217ab2063SBarry Smith   }
1453416022c9SBarry Smith   else c->diag          = 0;
14546c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
14556c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1456754ec7b1SSatish Balay   if (a->inode.size){
1457754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1458754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1459754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1460754ec7b1SSatish Balay   } else {
1461754ec7b1SSatish Balay     c->inode.size       = 0;
1462754ec7b1SSatish Balay     c->inode.node_count = 0;
1463754ec7b1SSatish Balay   }
1464416022c9SBarry Smith   c->nz                 = a->nz;
1465416022c9SBarry Smith   c->maxnz              = a->maxnz;
1466416022c9SBarry Smith   c->solve_work         = 0;
146776dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1468754ec7b1SSatish Balay 
1469416022c9SBarry Smith   *B = C;
147017ab2063SBarry Smith   return 0;
147117ab2063SBarry Smith }
147217ab2063SBarry Smith 
1473416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
147417ab2063SBarry Smith {
1475416022c9SBarry Smith   Mat_SeqAIJ   *a;
1476416022c9SBarry Smith   Mat          B;
147717699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
147817ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
147917ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
148017ab2063SBarry Smith 
148117699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
148217699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
148317ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1484416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
148548d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
148617ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
148717ab2063SBarry Smith 
148817ab2063SBarry Smith   /* read in row lengths */
14890452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1490416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
149117ab2063SBarry Smith 
149217ab2063SBarry Smith   /* create our matrix */
1493416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1494416022c9SBarry Smith   B = *A;
1495416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1496416022c9SBarry Smith   shift = a->indexshift;
149717ab2063SBarry Smith 
149817ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1499416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
150017ab2063SBarry Smith   if (shift) {
150117ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1502416022c9SBarry Smith       a->j[i] += 1;
150317ab2063SBarry Smith     }
150417ab2063SBarry Smith   }
150517ab2063SBarry Smith 
150617ab2063SBarry Smith   /* read in nonzero values */
1507416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
150817ab2063SBarry Smith 
150917ab2063SBarry Smith   /* set matrix "i" values */
1510416022c9SBarry Smith   a->i[0] = -shift;
151117ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1512416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1513416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
151417ab2063SBarry Smith   }
15150452661fSBarry Smith   PetscFree(rowlengths);
151617ab2063SBarry Smith 
1517416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1518416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
151917ab2063SBarry Smith   return 0;
152017ab2063SBarry Smith }
152117ab2063SBarry Smith 
152217ab2063SBarry Smith 
152317ab2063SBarry Smith 
15247264ac53SSatish Balay /*
15257264ac53SSatish Balay 
15267264ac53SSatish Balay   flg =0 if error;
15277264ac53SSatish Balay   flg =1 if both are equal;
15287264ac53SSatish Balay   */
15297264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg)
15307264ac53SSatish Balay {
15317264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
15327264ac53SSatish Balay 
15337264ac53SSatish Balay   if (B->type !=MATSEQAIJ)  SETERRQ(1,"MatEqual_SeqAIJ:Both matrices should be of type MATSEQAIJ");
15347264ac53SSatish Balay 
15357264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
15367264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
15377264ac53SSatish Balay       (a->indexshift != b->indexshift)) { *flg =0 ; return 0; }
15387264ac53SSatish Balay 
15397264ac53SSatish Balay   /* if the a->i are the same */
15407264ac53SSatish Balay   if(PetscMemcmp((char *)a->i, (char*)b->i, (a->n+1)*sizeof(int))) { *flg =0 ; return 0;}
15417264ac53SSatish Balay 
15427264ac53SSatish Balay   /* if a->j are the same */
15437264ac53SSatish Balay   if(PetscMemcmp((char *)a->j, (char*)b->j, (a->nz)*sizeof(int))) {
15447264ac53SSatish Balay     *flg =0 ; return 0;
15457264ac53SSatish Balay   }
15467264ac53SSatish Balay 
15477264ac53SSatish Balay   /* if a->j are the same */
15487264ac53SSatish Balay   if(PetscMemcmp((char *)a->a, (char*)b->a, (a->nz)*sizeof(Scalar))) {
15497264ac53SSatish Balay     *flg =0 ; return 0;
15507264ac53SSatish Balay   }
15517264ac53SSatish Balay   *flg =1 ;
15527264ac53SSatish Balay   return 0;
15537264ac53SSatish Balay 
15547264ac53SSatish Balay }
1555