xref: /petsc/src/mat/impls/aij/seq/aij.c (revision dfa27b742cf2ba09cbb22bc417bb67c000ed546d)
16d84be18SBarry Smith 
217ab2063SBarry Smith #ifndef lint
3*dfa27b74SSatish Balay static char vcid[] = "$Id: aij.c,v 1.166 1996/03/31 19:59:18 curfman Exp balay $";
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 
16bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,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;
21bcd2baecSBarry Smith   int        ierr, *ia, *ja,n,*idx,i,oshift,ishift;
2217ab2063SBarry Smith 
23a2744918SBarry Smith   /*
24a2744918SBarry Smith      this is tacky: In the future when we have written special factorization
25a2744918SBarry Smith      and solve routines for the identity permutation we should use a
26a2744918SBarry Smith      stride index set instead of the general one.
27a2744918SBarry Smith   */
28a2744918SBarry Smith   if (type  == ORDER_NATURAL) {
29a2744918SBarry Smith     n = a->n;
30a2744918SBarry Smith     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
31a2744918SBarry Smith     for ( i=0; i<n; i++ ) idx[i] = i;
32a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
33a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
34a2744918SBarry Smith     PetscFree(idx);
35a2744918SBarry Smith     ISSetPermutation(*rperm);
36a2744918SBarry Smith     ISSetPermutation(*cperm);
37a2744918SBarry Smith     ISSetIdentity(*rperm);
38a2744918SBarry Smith     ISSetIdentity(*cperm);
39a2744918SBarry Smith     return 0;
40a2744918SBarry Smith   }
41a2744918SBarry Smith 
42bcd2baecSBarry Smith   MatReorderingRegisterAll();
43bcd2baecSBarry Smith   ishift = a->indexshift;
44bcd2baecSBarry Smith   oshift = -MatReorderingIndexShift[(int)type];
45bcd2baecSBarry Smith   if (MatReorderingRequiresSymmetric[(int)type]) {
46bcd2baecSBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja);
47bcd2baecSBarry Smith     CHKERRQ(ierr);
48416022c9SBarry Smith     ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
490452661fSBarry Smith     PetscFree(ia); PetscFree(ja);
50bcd2baecSBarry Smith   } else {
51bcd2baecSBarry Smith     if (ishift == oshift) {
52bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
53bcd2baecSBarry Smith     }
54bcd2baecSBarry Smith     else if (ishift == -1) {
55bcd2baecSBarry Smith       /* temporarily subtract 1 from i and j indices */
56bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
57bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
58bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
59bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
60bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
61bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
62bcd2baecSBarry Smith     } else {
63bcd2baecSBarry Smith       /* temporarily add 1 to i and j indices */
64bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
65bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
66bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
67bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
68bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
69bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
70bcd2baecSBarry Smith     }
71bcd2baecSBarry Smith   }
7217ab2063SBarry Smith   return 0;
7317ab2063SBarry Smith }
7417ab2063SBarry Smith 
75227d817aSBarry Smith #define CHUNKSIZE   15
7617ab2063SBarry Smith 
7717ab2063SBarry Smith /* This version has row oriented v  */
78416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
7917ab2063SBarry Smith {
80416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
81416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
824b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
83d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
84416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
8517ab2063SBarry Smith 
8617ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
87416022c9SBarry Smith     row  = im[k];
8817ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
89416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
9017ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
9117ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
92416022c9SBarry Smith     low = 0;
9317ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
94416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
95416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
964b0e389bSBarry Smith       col = in[l] - shift;
974b0e389bSBarry Smith       if (roworiented) {
984b0e389bSBarry Smith         value = *v++;
994b0e389bSBarry Smith       }
1004b0e389bSBarry Smith       else {
1014b0e389bSBarry Smith         value = v[k + l*m];
1024b0e389bSBarry Smith       }
103416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
104416022c9SBarry Smith       while (high-low > 5) {
105416022c9SBarry Smith         t = (low+high)/2;
106416022c9SBarry Smith         if (rp[t] > col) high = t;
107416022c9SBarry Smith         else             low  = t;
10817ab2063SBarry Smith       }
109416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
11017ab2063SBarry Smith         if (rp[i] > col) break;
11117ab2063SBarry Smith         if (rp[i] == col) {
112416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
11317ab2063SBarry Smith           else                  ap[i] = value;
11417ab2063SBarry Smith           goto noinsert;
11517ab2063SBarry Smith         }
11617ab2063SBarry Smith       }
11717ab2063SBarry Smith       if (nonew) goto noinsert;
11817ab2063SBarry Smith       if (nrow >= rmax) {
11917ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
120416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
12117ab2063SBarry Smith         Scalar *new_a;
12217ab2063SBarry Smith 
12317ab2063SBarry Smith         /* malloc new storage space */
124416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1250452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
12617ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
12717ab2063SBarry Smith         new_i   = new_j + new_nz;
12817ab2063SBarry Smith 
12917ab2063SBarry Smith         /* copy over old data into new slots */
13017ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
131416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
132416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
133416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
134416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
13517ab2063SBarry Smith                                                            len*sizeof(int));
136416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
137416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
13817ab2063SBarry Smith                                                            len*sizeof(Scalar));
13917ab2063SBarry Smith         /* free up old matrix storage */
1400452661fSBarry Smith         PetscFree(a->a);
1410452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
142416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
143416022c9SBarry Smith         a->singlemalloc = 1;
14417ab2063SBarry Smith 
14517ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
146416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
147416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
148416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
149b810aeb4SBarry Smith         a->reallocs++;
15017ab2063SBarry Smith       }
151416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
152416022c9SBarry Smith       /* shift up all the later entries in this row */
153416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
15417ab2063SBarry Smith         rp[ii+1] = rp[ii];
15517ab2063SBarry Smith         ap[ii+1] = ap[ii];
15617ab2063SBarry Smith       }
15717ab2063SBarry Smith       rp[i] = col;
15817ab2063SBarry Smith       ap[i] = value;
15917ab2063SBarry Smith       noinsert:;
160416022c9SBarry Smith       low = i + 1;
16117ab2063SBarry Smith     }
16217ab2063SBarry Smith     ailen[row] = nrow;
16317ab2063SBarry Smith   }
16417ab2063SBarry Smith   return 0;
16517ab2063SBarry Smith }
16617ab2063SBarry Smith 
1677eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1687eb43aa7SLois Curfman McInnes {
1697eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
170b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1717eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1727eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
1737eb43aa7SLois Curfman McInnes 
1747eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
1757eb43aa7SLois Curfman McInnes     row  = im[k];
1767eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
1777eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
1787eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
1797eb43aa7SLois Curfman McInnes     nrow = ailen[row];
1807eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
1817eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
1827eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
1837eb43aa7SLois Curfman McInnes       col = in[l] - shift;
1847eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
1857eb43aa7SLois Curfman McInnes       while (high-low > 5) {
1867eb43aa7SLois Curfman McInnes         t = (low+high)/2;
1877eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
1887eb43aa7SLois Curfman McInnes         else             low  = t;
1897eb43aa7SLois Curfman McInnes       }
1907eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
1917eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
1927eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
193b49de8d1SLois Curfman McInnes           *v++ = ap[i];
1947eb43aa7SLois Curfman McInnes           goto finished;
1957eb43aa7SLois Curfman McInnes         }
1967eb43aa7SLois Curfman McInnes       }
197b49de8d1SLois Curfman McInnes       *v++ = zero;
1987eb43aa7SLois Curfman McInnes       finished:;
1997eb43aa7SLois Curfman McInnes     }
2007eb43aa7SLois Curfman McInnes   }
2017eb43aa7SLois Curfman McInnes   return 0;
2027eb43aa7SLois Curfman McInnes }
2037eb43aa7SLois Curfman McInnes 
20417ab2063SBarry Smith #include "draw.h"
20517ab2063SBarry Smith #include "pinclude/pviewer.h"
20677c4ece6SBarry Smith #include "sys.h"
20717ab2063SBarry Smith 
208416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
20917ab2063SBarry Smith {
210416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
211416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
21217ab2063SBarry Smith 
21390ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2140452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
215416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
216416022c9SBarry Smith   col_lens[1] = a->m;
217416022c9SBarry Smith   col_lens[2] = a->n;
218416022c9SBarry Smith   col_lens[3] = a->nz;
219416022c9SBarry Smith 
220416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
221416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
222416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
22317ab2063SBarry Smith   }
22477c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
2250452661fSBarry Smith   PetscFree(col_lens);
226416022c9SBarry Smith 
227416022c9SBarry Smith   /* store column indices (zero start index) */
228416022c9SBarry Smith   if (a->indexshift) {
229416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
23017ab2063SBarry Smith   }
23177c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr);
232416022c9SBarry Smith   if (a->indexshift) {
233416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
23417ab2063SBarry Smith   }
235416022c9SBarry Smith 
236416022c9SBarry Smith   /* store nonzero values */
23777c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
23817ab2063SBarry Smith   return 0;
23917ab2063SBarry Smith }
240416022c9SBarry Smith 
241416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
242416022c9SBarry Smith {
243416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
24495e01e2fSLois Curfman McInnes   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg;
24517ab2063SBarry Smith   FILE        *fd;
24617ab2063SBarry Smith   char        *outputname;
24717ab2063SBarry Smith 
24890ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
249416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
25090ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
25190ace30eSBarry Smith   if (format == ASCII_FORMAT_INFO) {
25295e01e2fSLois Curfman McInnes     return 0;
25395e01e2fSLois Curfman McInnes   }
25490ace30eSBarry Smith   else if (format == ASCII_FORMAT_INFO_DETAILED) {
25595e01e2fSLois Curfman McInnes     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr);
25695e01e2fSLois Curfman McInnes     if (flg) fprintf(fd,"  not using I-node routines\n");
25795e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
25895e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
25917ab2063SBarry Smith   }
26090ace30eSBarry Smith   else if (format == ASCII_FORMAT_MATLAB) {
26117ab2063SBarry Smith     int nz, nzalloc, mem;
262416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
263416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
26417ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
26517ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
26617ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
26717ab2063SBarry Smith 
26817ab2063SBarry Smith     for (i=0; i<m; i++) {
269416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
27017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
2717a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]),
272416022c9SBarry Smith                    imag(a->a[j]));
27317ab2063SBarry Smith #else
2747a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
27517ab2063SBarry Smith #endif
27617ab2063SBarry Smith       }
27717ab2063SBarry Smith     }
27817ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
27917ab2063SBarry Smith   }
28017ab2063SBarry Smith   else {
28117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
28217ab2063SBarry Smith       fprintf(fd,"row %d:",i);
283416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
28417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
285416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
286416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
28717ab2063SBarry Smith         }
28817ab2063SBarry Smith         else {
289416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
29017ab2063SBarry Smith         }
29117ab2063SBarry Smith #else
292416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
29317ab2063SBarry Smith #endif
29417ab2063SBarry Smith       }
29517ab2063SBarry Smith       fprintf(fd,"\n");
29617ab2063SBarry Smith     }
29717ab2063SBarry Smith   }
29817ab2063SBarry Smith   fflush(fd);
299416022c9SBarry Smith   return 0;
300416022c9SBarry Smith }
301416022c9SBarry Smith 
302416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
303416022c9SBarry Smith {
304416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
305cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
306cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
307bcd2baecSBarry Smith   Draw        draw;
308cddf8d76SBarry Smith   DrawButton  button;
30919bcc07fSBarry Smith   PetscTruth  isnull;
310cddf8d76SBarry Smith 
311bcd2baecSBarry Smith   ViewerDrawGetDraw(viewer,&draw);
31219bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
31319bcc07fSBarry Smith 
314416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
315416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
316416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
317416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
318cddf8d76SBarry Smith   color = DRAW_BLUE;
319416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
320cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
321416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
322cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
323cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
324cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
325cddf8d76SBarry Smith #else
326cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
327cddf8d76SBarry Smith #endif
328cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
329cddf8d76SBarry Smith     }
330cddf8d76SBarry Smith   }
331cddf8d76SBarry Smith   color = DRAW_CYAN;
332cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
333cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
334cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
335cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
336cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
337cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
338cddf8d76SBarry Smith     }
339cddf8d76SBarry Smith   }
340cddf8d76SBarry Smith   color = DRAW_RED;
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);
351416022c9SBarry Smith     }
352416022c9SBarry Smith   }
353416022c9SBarry Smith   DrawFlush(draw);
354cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
355cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
356cddf8d76SBarry Smith 
357cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
358cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
359cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
360cddf8d76SBarry Smith     DrawClear(draw);
361cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
362cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
363cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
364cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
365cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
366cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
367cddf8d76SBarry Smith     w *= scale; h *= scale;
368cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
369cddf8d76SBarry Smith     color = DRAW_BLUE;
370cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
371cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
372cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
373cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
374cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
375cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
376cddf8d76SBarry Smith #else
377cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
378cddf8d76SBarry Smith #endif
379cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
380cddf8d76SBarry Smith       }
381cddf8d76SBarry Smith     }
382cddf8d76SBarry Smith     color = DRAW_CYAN;
383cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
384cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
385cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
386cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
387cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
388cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
389cddf8d76SBarry Smith       }
390cddf8d76SBarry Smith     }
391cddf8d76SBarry Smith     color = DRAW_RED;
392cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
393cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
394cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
395cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
396cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
397cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
398cddf8d76SBarry Smith #else
399cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
400cddf8d76SBarry Smith #endif
401cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
402cddf8d76SBarry Smith       }
403cddf8d76SBarry Smith     }
404cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
405cddf8d76SBarry Smith   }
406416022c9SBarry Smith   return 0;
407416022c9SBarry Smith }
408416022c9SBarry Smith 
409416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
410416022c9SBarry Smith {
411416022c9SBarry Smith   Mat         A = (Mat) obj;
412416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
413bcd2baecSBarry Smith   ViewerType  vtype;
414bcd2baecSBarry Smith   int         ierr;
415416022c9SBarry Smith 
416416022c9SBarry Smith   if (!viewer) {
417bcd2baecSBarry Smith     viewer = STDOUT_VIEWER_SELF;
418416022c9SBarry Smith   }
419bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
420bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
421416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
422416022c9SBarry Smith   }
423bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
424416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
425416022c9SBarry Smith   }
426bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
427416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
428416022c9SBarry Smith   }
429bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
430bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
43117ab2063SBarry Smith   }
43217ab2063SBarry Smith   return 0;
43317ab2063SBarry Smith }
43419bcc07fSBarry Smith 
435c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
436416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
43717ab2063SBarry Smith {
438416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
43941c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
440416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
441416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
44217ab2063SBarry Smith 
44317ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
44417ab2063SBarry Smith 
44517ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
446416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
44717ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
44817ab2063SBarry Smith     if (fshift) {
449416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
45017ab2063SBarry Smith       N = ailen[i];
45117ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
45217ab2063SBarry Smith         ip[j-fshift] = ip[j];
45317ab2063SBarry Smith         ap[j-fshift] = ap[j];
45417ab2063SBarry Smith       }
45517ab2063SBarry Smith     }
45617ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
45717ab2063SBarry Smith   }
45817ab2063SBarry Smith   if (m) {
45917ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
46017ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
46117ab2063SBarry Smith   }
46217ab2063SBarry Smith   /* reset ilen and imax for each row */
46317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
46417ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
46517ab2063SBarry Smith   }
466416022c9SBarry Smith   a->nz = ai[m] + shift;
46717ab2063SBarry Smith 
46817ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
469416022c9SBarry Smith   if (fshift && a->diag) {
4700452661fSBarry Smith     PetscFree(a->diag);
471416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
472416022c9SBarry Smith     a->diag = 0;
47317ab2063SBarry Smith   }
47494a424c1SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n",
475b810aeb4SBarry Smith            fshift,a->nz,m);
47694a424c1SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n",
477b810aeb4SBarry Smith            a->reallocs);
47876dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
47941c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
48017ab2063SBarry Smith   return 0;
48117ab2063SBarry Smith }
48217ab2063SBarry Smith 
483416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
48417ab2063SBarry Smith {
485416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
486cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
48717ab2063SBarry Smith   return 0;
48817ab2063SBarry Smith }
489416022c9SBarry Smith 
49017ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
49117ab2063SBarry Smith {
492416022c9SBarry Smith   Mat        A  = (Mat) obj;
493416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
494d5d45c9bSBarry Smith 
49517ab2063SBarry Smith #if defined(PETSC_LOG)
496416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
49717ab2063SBarry Smith #endif
4980452661fSBarry Smith   PetscFree(a->a);
4990452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
5000452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
5010452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
5020452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
5030452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
50476dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
5050452661fSBarry Smith   PetscFree(a);
506f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
507f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
50817ab2063SBarry Smith   return 0;
50917ab2063SBarry Smith }
51017ab2063SBarry Smith 
511416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
51217ab2063SBarry Smith {
51317ab2063SBarry Smith   return 0;
51417ab2063SBarry Smith }
51517ab2063SBarry Smith 
516416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
51717ab2063SBarry Smith {
518416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
519416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
5204b0e389bSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
521416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
522416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
523416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
524e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
525e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
526e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
527e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
52894a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
529df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
5304d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
5316c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
5326c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
5336c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
5346c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
5356c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
536e2f28af5SBarry Smith   else
5374d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
53817ab2063SBarry Smith   return 0;
53917ab2063SBarry Smith }
54017ab2063SBarry Smith 
541416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
54217ab2063SBarry Smith {
543416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
544416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
54517ab2063SBarry Smith   Scalar     *x, zero = 0.0;
54617ab2063SBarry Smith 
54717ab2063SBarry Smith   VecSet(&zero,v);
54817ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
549416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
550416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
551416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
552416022c9SBarry Smith       if (a->j[j]+shift == i) {
553416022c9SBarry Smith         x[i] = a->a[j];
55417ab2063SBarry Smith         break;
55517ab2063SBarry Smith       }
55617ab2063SBarry Smith     }
55717ab2063SBarry Smith   }
55817ab2063SBarry Smith   return 0;
55917ab2063SBarry Smith }
56017ab2063SBarry Smith 
56117ab2063SBarry Smith /* -------------------------------------------------------*/
56217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
56317ab2063SBarry Smith /* -------------------------------------------------------*/
564416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
56517ab2063SBarry Smith {
566416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
56717ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
568416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
56917ab2063SBarry Smith 
57017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
571cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
57217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
57317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
574416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
575416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
576416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
57717ab2063SBarry Smith     alpha = x[i];
57817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
57917ab2063SBarry Smith   }
580416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
58117ab2063SBarry Smith   return 0;
58217ab2063SBarry Smith }
583d5d45c9bSBarry Smith 
584416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
58517ab2063SBarry Smith {
586416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
58717ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
588416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
58917ab2063SBarry Smith 
59017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
59117ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
59217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
59317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
594416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
595416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
596416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
59717ab2063SBarry Smith     alpha = x[i];
59817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
59917ab2063SBarry Smith   }
60017ab2063SBarry Smith   return 0;
60117ab2063SBarry Smith }
60217ab2063SBarry Smith 
603416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
60417ab2063SBarry Smith {
605416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
60617ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
607416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
60817ab2063SBarry Smith 
60917ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
61017ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
611416022c9SBarry Smith   idx  = a->j;
612416022c9SBarry Smith   v    = a->a;
613416022c9SBarry Smith   ii   = a->i;
61417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
615416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
61617ab2063SBarry Smith     sum  = 0.0;
61717ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
61817ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
619416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
62017ab2063SBarry Smith     y[i] = sum;
62117ab2063SBarry Smith   }
622416022c9SBarry Smith   PLogFlops(2*a->nz - m);
62317ab2063SBarry Smith   return 0;
62417ab2063SBarry Smith }
62517ab2063SBarry Smith 
626416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
62717ab2063SBarry Smith {
628416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
62917ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
630cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
63117ab2063SBarry Smith 
63217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
63317ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
634cddf8d76SBarry Smith   idx  = a->j;
635cddf8d76SBarry Smith   v    = a->a;
636cddf8d76SBarry Smith   ii   = a->i;
63717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
638cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
63917ab2063SBarry Smith     sum  = y[i];
640cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
641cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
64217ab2063SBarry Smith     z[i] = sum;
64317ab2063SBarry Smith   }
644416022c9SBarry Smith   PLogFlops(2*a->nz);
64517ab2063SBarry Smith   return 0;
64617ab2063SBarry Smith }
64717ab2063SBarry Smith 
64817ab2063SBarry Smith /*
64917ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
65017ab2063SBarry Smith */
65117ab2063SBarry Smith 
652416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
65317ab2063SBarry Smith {
654416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
655416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
65617ab2063SBarry Smith 
6570452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
658416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
659416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
660416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
661416022c9SBarry Smith       if (a->j[j]+shift == i) {
66217ab2063SBarry Smith         diag[i] = j - shift;
66317ab2063SBarry Smith         break;
66417ab2063SBarry Smith       }
66517ab2063SBarry Smith     }
66617ab2063SBarry Smith   }
667416022c9SBarry Smith   a->diag = diag;
66817ab2063SBarry Smith   return 0;
66917ab2063SBarry Smith }
67017ab2063SBarry Smith 
671416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
67217ab2063SBarry Smith                            double fshift,int its,Vec xx)
67317ab2063SBarry Smith {
674416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
675416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
676d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
67717ab2063SBarry Smith 
67817ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
679416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
680416022c9SBarry Smith   diag = a->diag;
681416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
68217ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
68317ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
68417ab2063SBarry Smith     bs = b + shift;
68517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
686416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
687416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
688416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
689416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
69017ab2063SBarry Smith         sum  = b[i]*d/omega;
69117ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
69217ab2063SBarry Smith         x[i] = sum;
69317ab2063SBarry Smith     }
69417ab2063SBarry Smith     return 0;
69517ab2063SBarry Smith   }
69617ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
69717ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
69817ab2063SBarry Smith   }
699416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
70017ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
70117ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
70217ab2063SBarry Smith 
70317ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
70417ab2063SBarry Smith 
70517ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
70617ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
70717ab2063SBarry Smith     is the relaxation factor.
70817ab2063SBarry Smith     */
7090452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
71017ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
71117ab2063SBarry Smith 
71217ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
71317ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
714416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
715416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
716416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
717416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
71817ab2063SBarry Smith       sum  = b[i];
71917ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
72017ab2063SBarry Smith       x[i] = omega*(sum/d);
72117ab2063SBarry Smith     }
72217ab2063SBarry Smith 
72317ab2063SBarry Smith     /*  t = b - (2*E - D)x */
724416022c9SBarry Smith     v = a->a;
72517ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
72617ab2063SBarry Smith 
72717ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
728416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
729416022c9SBarry Smith     diag = a->diag;
73017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
731416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
732416022c9SBarry Smith       n    = diag[i] - a->i[i];
733416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
734416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
73517ab2063SBarry Smith       sum  = t[i];
73617ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
73717ab2063SBarry Smith       t[i] = omega*(sum/d);
73817ab2063SBarry Smith     }
73917ab2063SBarry Smith 
74017ab2063SBarry Smith     /*  x = x + t */
74117ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7420452661fSBarry Smith     PetscFree(t);
74317ab2063SBarry Smith     return 0;
74417ab2063SBarry Smith   }
74517ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
74617ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
74717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
748416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
749416022c9SBarry Smith         n    = diag[i] - a->i[i];
750416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
751416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
75217ab2063SBarry Smith         sum  = b[i];
75317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
75417ab2063SBarry Smith         x[i] = omega*(sum/d);
75517ab2063SBarry Smith       }
75617ab2063SBarry Smith       xb = x;
75717ab2063SBarry Smith     }
75817ab2063SBarry Smith     else xb = b;
75917ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
76017ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
76117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
762416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
76317ab2063SBarry Smith       }
76417ab2063SBarry Smith     }
76517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
76617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
767416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
768416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
769416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
770416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
77117ab2063SBarry Smith         sum  = xb[i];
77217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
77317ab2063SBarry Smith         x[i] = omega*(sum/d);
77417ab2063SBarry Smith       }
77517ab2063SBarry Smith     }
77617ab2063SBarry Smith     its--;
77717ab2063SBarry Smith   }
77817ab2063SBarry Smith   while (its--) {
77917ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
78017ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
781416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
782416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
783416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
784416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
78517ab2063SBarry Smith         sum  = b[i];
78617ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
78717ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
78817ab2063SBarry Smith       }
78917ab2063SBarry Smith     }
79017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
79117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
792416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
793416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
794416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
795416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
79617ab2063SBarry Smith         sum  = b[i];
79717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
79817ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
79917ab2063SBarry Smith       }
80017ab2063SBarry Smith     }
80117ab2063SBarry Smith   }
80217ab2063SBarry Smith   return 0;
80317ab2063SBarry Smith }
80417ab2063SBarry Smith 
805d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
80617ab2063SBarry Smith {
807416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
808bcd2baecSBarry Smith   if (nz)      *nz      = a->nz;
809bcd2baecSBarry Smith   if (nzalloc) *nzalloc = a->maxnz;
810bcd2baecSBarry Smith   if (mem)     *mem     = (int)A->mem;
81117ab2063SBarry Smith   return 0;
81217ab2063SBarry Smith }
81317ab2063SBarry Smith 
81417ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
81517ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
81617ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
81717ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
81817ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
81917ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
82017ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
82117ab2063SBarry Smith 
82217ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
82317ab2063SBarry Smith {
824416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
825416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
82617ab2063SBarry Smith 
82777c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
82817ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
82917ab2063SBarry Smith   if (diag) {
83017ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
831416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
832416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
833416022c9SBarry Smith         a->ilen[rows[i]] = 1;
834416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
835416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
83617ab2063SBarry Smith       }
83717ab2063SBarry Smith       else {
83817ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
83917ab2063SBarry Smith         CHKERRQ(ierr);
84017ab2063SBarry Smith       }
84117ab2063SBarry Smith     }
84217ab2063SBarry Smith   }
84317ab2063SBarry Smith   else {
84417ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
845416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
846416022c9SBarry Smith       a->ilen[rows[i]] = 0;
84717ab2063SBarry Smith     }
84817ab2063SBarry Smith   }
84917ab2063SBarry Smith   ISRestoreIndices(is,&rows);
85017ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
85117ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
85217ab2063SBarry Smith   return 0;
85317ab2063SBarry Smith }
85417ab2063SBarry Smith 
855416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
85617ab2063SBarry Smith {
857416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
858416022c9SBarry Smith   *m = a->m; *n = a->n;
85917ab2063SBarry Smith   return 0;
86017ab2063SBarry Smith }
86117ab2063SBarry Smith 
862416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
86317ab2063SBarry Smith {
864416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
865416022c9SBarry Smith   *m = 0; *n = a->m;
86617ab2063SBarry Smith   return 0;
86717ab2063SBarry Smith }
8684e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
86917ab2063SBarry Smith {
870416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
871c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
87217ab2063SBarry Smith 
873416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
87417ab2063SBarry Smith 
875416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
876416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
87717ab2063SBarry Smith   if (idx) {
878416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
8794e093b46SBarry Smith     if (*nz && shift) {
8800452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
88117ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
8824e093b46SBarry Smith     } else if (*nz) {
8834e093b46SBarry Smith       *idx = itmp;
88417ab2063SBarry Smith     }
88517ab2063SBarry Smith     else *idx = 0;
88617ab2063SBarry Smith   }
88717ab2063SBarry Smith   return 0;
88817ab2063SBarry Smith }
88917ab2063SBarry Smith 
8904e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
89117ab2063SBarry Smith {
8924e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
8934e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
89417ab2063SBarry Smith   return 0;
89517ab2063SBarry Smith }
89617ab2063SBarry Smith 
897cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
89817ab2063SBarry Smith {
899416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
900416022c9SBarry Smith   Scalar     *v = a->a;
90117ab2063SBarry Smith   double     sum = 0.0;
902416022c9SBarry Smith   int        i, j,shift = a->indexshift;
90317ab2063SBarry Smith 
90417ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
905416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
90617ab2063SBarry Smith #if defined(PETSC_COMPLEX)
90717ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
90817ab2063SBarry Smith #else
90917ab2063SBarry Smith       sum += (*v)*(*v); v++;
91017ab2063SBarry Smith #endif
91117ab2063SBarry Smith     }
91217ab2063SBarry Smith     *norm = sqrt(sum);
91317ab2063SBarry Smith   }
91417ab2063SBarry Smith   else if (type == NORM_1) {
91517ab2063SBarry Smith     double *tmp;
916416022c9SBarry Smith     int    *jj = a->j;
9170452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
918cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
91917ab2063SBarry Smith     *norm = 0.0;
920416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
921a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
92217ab2063SBarry Smith     }
923416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
92417ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
92517ab2063SBarry Smith     }
9260452661fSBarry Smith     PetscFree(tmp);
92717ab2063SBarry Smith   }
92817ab2063SBarry Smith   else if (type == NORM_INFINITY) {
92917ab2063SBarry Smith     *norm = 0.0;
930416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
931416022c9SBarry Smith       v = a->a + a->i[j] + shift;
93217ab2063SBarry Smith       sum = 0.0;
933416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
934cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
93517ab2063SBarry Smith       }
93617ab2063SBarry Smith       if (sum > *norm) *norm = sum;
93717ab2063SBarry Smith     }
93817ab2063SBarry Smith   }
93917ab2063SBarry Smith   else {
94048d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
94117ab2063SBarry Smith   }
94217ab2063SBarry Smith   return 0;
94317ab2063SBarry Smith }
94417ab2063SBarry Smith 
945416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
94617ab2063SBarry Smith {
947416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
948416022c9SBarry Smith   Mat        C;
949416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
950416022c9SBarry Smith   int        shift = a->indexshift;
951d5d45c9bSBarry Smith   Scalar     *array = a->a;
95217ab2063SBarry Smith 
9533638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
954*dfa27b74SSatish Balay     SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place");
9550452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
956cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
95717ab2063SBarry Smith   if (shift) {
95817ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
95917ab2063SBarry Smith   }
96017ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
961416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
9620452661fSBarry Smith   PetscFree(col);
96317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
96417ab2063SBarry Smith     len = ai[i+1]-ai[i];
965416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
96617ab2063SBarry Smith     array += len; aj += len;
96717ab2063SBarry Smith   }
96817ab2063SBarry Smith   if (shift) {
96917ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
97017ab2063SBarry Smith   }
97117ab2063SBarry Smith 
972416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
973416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
97417ab2063SBarry Smith 
9753638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
976416022c9SBarry Smith     *B = C;
97717ab2063SBarry Smith   } else {
978416022c9SBarry Smith     /* This isn't really an in-place transpose */
9790452661fSBarry Smith     PetscFree(a->a);
9800452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
9810452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
9820452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
9830452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
9840452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
9850452661fSBarry Smith     PetscFree(a);
986416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
9870452661fSBarry Smith     PetscHeaderDestroy(C);
98817ab2063SBarry Smith   }
98917ab2063SBarry Smith   return 0;
99017ab2063SBarry Smith }
99117ab2063SBarry Smith 
992f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
99317ab2063SBarry Smith {
994416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
99517ab2063SBarry Smith   Scalar     *l,*r,x,*v;
996d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
99717ab2063SBarry Smith 
99817ab2063SBarry Smith   if (ll) {
99917ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
1000f0b747eeSBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
1001416022c9SBarry Smith     v = a->a;
100217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
100317ab2063SBarry Smith       x = l[i];
1004416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
100517ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
100617ab2063SBarry Smith     }
100717ab2063SBarry Smith   }
100817ab2063SBarry Smith   if (rr) {
100917ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
1010f0b747eeSBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
1011416022c9SBarry Smith     v = a->a; jj = a->j;
101217ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
101317ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
101417ab2063SBarry Smith     }
101517ab2063SBarry Smith   }
101617ab2063SBarry Smith   return 0;
101717ab2063SBarry Smith }
101817ab2063SBarry Smith 
1019cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
102017ab2063SBarry Smith {
1021db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
102202834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
102399141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1024a2744918SBarry Smith   register int sum,lensi;
102599141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
102699141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
102799141d43SSatish Balay   Scalar       *a_new,*mat_a;
1028416022c9SBarry Smith   Mat          C;
102917ab2063SBarry Smith 
103099141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
103199141d43SSatish Balay   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IS is not sorted");
103299141d43SSatish Balay 
103317ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
103417ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
103517ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
103617ab2063SBarry Smith 
10377264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
103802834360SBarry Smith     /* special case of contiguous rows */
103957faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
104002834360SBarry Smith     starts = lens + ncols;
104102834360SBarry Smith     /* loop over new rows determining lens and starting points */
104202834360SBarry Smith     for (i=0; i<nrows; i++) {
1043a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1044a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
104502834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1046d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
104702834360SBarry Smith           starts[i] = k;
104802834360SBarry Smith           break;
104902834360SBarry Smith 	}
105002834360SBarry Smith       }
1051a2744918SBarry Smith       sum = 0;
105202834360SBarry Smith       while (k < kend) {
1053d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1054a2744918SBarry Smith         sum++;
105502834360SBarry Smith       }
1056a2744918SBarry Smith       lens[i] = sum;
105702834360SBarry Smith     }
105802834360SBarry Smith     /* create submatrix */
1059cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
106008480c60SBarry Smith       int n_cols,n_rows;
106108480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
106208480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1063d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
106408480c60SBarry Smith       C = *B;
106508480c60SBarry Smith     }
106608480c60SBarry Smith     else {
106702834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
106808480c60SBarry Smith     }
1069db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1070db02288aSLois Curfman McInnes 
107102834360SBarry Smith     /* loop over rows inserting into submatrix */
1072db02288aSLois Curfman McInnes     a_new    = c->a;
1073db02288aSLois Curfman McInnes     j_new    = c->j;
1074db02288aSLois Curfman McInnes     i_new    = c->i;
1075db02288aSLois Curfman McInnes     i_new[0] = -shift;
107602834360SBarry Smith     for (i=0; i<nrows; i++) {
1077a2744918SBarry Smith       ii    = starts[i];
1078a2744918SBarry Smith       lensi = lens[i];
1079a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1080a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
108102834360SBarry Smith       }
1082a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1083a2744918SBarry Smith       a_new      += lensi;
1084a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1085a2744918SBarry Smith       c->ilen[i]  = lensi;
108602834360SBarry Smith     }
10870452661fSBarry Smith     PetscFree(lens);
108802834360SBarry Smith   }
108902834360SBarry Smith   else {
109002834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10910452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
109202834360SBarry Smith     ssmap = smap + shift;
109399141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1094cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
109517ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
109602834360SBarry Smith     /* determine lens of each row */
109702834360SBarry Smith     for (i=0; i<nrows; i++) {
1098d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
109902834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
110002834360SBarry Smith       lens[i] = 0;
110102834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1102d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
110302834360SBarry Smith           lens[i]++;
110402834360SBarry Smith         }
110502834360SBarry Smith       }
110602834360SBarry Smith     }
110717ab2063SBarry Smith     /* Create and fill new matrix */
1108a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
110999141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
111099141d43SSatish Balay 
111199141d43SSatish Balay       if (c->m  != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
111299141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
111399141d43SSatish Balay         SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros");
111499141d43SSatish Balay       }
111599141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
111608480c60SBarry Smith       C = *B;
111708480c60SBarry Smith     }
111808480c60SBarry Smith     else {
111902834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
112008480c60SBarry Smith     }
112199141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
112217ab2063SBarry Smith     for (i=0; i<nrows; i++) {
112399141d43SSatish Balay       row    = irow[i];
112417ab2063SBarry Smith       nznew  = 0;
112599141d43SSatish Balay       kstart = ai[row]+shift;
112699141d43SSatish Balay       kend   = kstart + a->ilen[row];
112799141d43SSatish Balay       mat_i  = c->i[i]+shift;
112899141d43SSatish Balay       mat_j  = c->j + mat_i;
112999141d43SSatish Balay       mat_a  = c->a + mat_i;
113099141d43SSatish Balay       mat_ilen = c->ilen + i;
113117ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
113299141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
113399141d43SSatish Balay           *mat_j++ = tcol - (!shift);
113499141d43SSatish Balay           *mat_a++ = a->a[k];
113599141d43SSatish Balay           (*mat_ilen)++;
113699141d43SSatish Balay 
113717ab2063SBarry Smith         }
113817ab2063SBarry Smith       }
113917ab2063SBarry Smith     }
114002834360SBarry Smith     /* Free work space */
114102834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
114299141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
114302834360SBarry Smith   }
1144416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1145416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
114617ab2063SBarry Smith 
114717ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1148416022c9SBarry Smith   *B = C;
114917ab2063SBarry Smith   return 0;
115017ab2063SBarry Smith }
115117ab2063SBarry Smith 
1152a871dcd8SBarry Smith /*
115363b91edcSBarry Smith      note: This can only work for identity for row and col. It would
115463b91edcSBarry Smith    be good to check this and otherwise generate an error.
1155a871dcd8SBarry Smith */
115663b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1157a871dcd8SBarry Smith {
115863b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
115908480c60SBarry Smith   int        ierr;
116063b91edcSBarry Smith   Mat        outA;
116163b91edcSBarry Smith 
1162a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1163a871dcd8SBarry Smith 
116463b91edcSBarry Smith   outA          = inA;
116563b91edcSBarry Smith   inA->factor   = FACTOR_LU;
116663b91edcSBarry Smith   a->row        = row;
116763b91edcSBarry Smith   a->col        = col;
116863b91edcSBarry Smith 
11690452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
117063b91edcSBarry Smith 
117108480c60SBarry Smith   if (!a->diag) {
117208480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
117363b91edcSBarry Smith   }
117463b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1175a871dcd8SBarry Smith   return 0;
1176a871dcd8SBarry Smith }
1177a871dcd8SBarry Smith 
1178f0b747eeSBarry Smith #include "pinclude/plapack.h"
1179f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1180f0b747eeSBarry Smith {
1181f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1182f0b747eeSBarry Smith   int        one = 1;
1183f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1184f0b747eeSBarry Smith   PLogFlops(a->nz);
1185f0b747eeSBarry Smith   return 0;
1186f0b747eeSBarry Smith }
1187f0b747eeSBarry Smith 
1188cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1189cddf8d76SBarry Smith                                     Mat **B)
1190cddf8d76SBarry Smith {
1191cddf8d76SBarry Smith   int ierr,i;
1192cddf8d76SBarry Smith 
1193cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
11940452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1195cddf8d76SBarry Smith   }
1196cddf8d76SBarry Smith 
1197cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1198cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1199cddf8d76SBarry Smith   }
1200cddf8d76SBarry Smith   return 0;
1201cddf8d76SBarry Smith }
1202cddf8d76SBarry Smith 
1203e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
12044dcbc457SBarry Smith {
1205e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
120606763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
12078a047759SSatish Balay   int        start, end, *ai, *aj;
120806763907SSatish Balay   char       *table;
12098a047759SSatish Balay   shift = a->indexshift;
1210e4d965acSSatish Balay   m     = a->m;
1211e4d965acSSatish Balay   ai    = a->i;
12128a047759SSatish Balay   aj    = a->j+shift;
12138a047759SSatish Balay 
1214e4d965acSSatish Balay   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
121506763907SSatish Balay 
121606763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
121706763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
121806763907SSatish Balay 
1219e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1220e4d965acSSatish Balay     /* Initialise the two local arrays */
1221e4d965acSSatish Balay     isz  = 0;
122206763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1223e4d965acSSatish Balay 
1224e4d965acSSatish Balay                 /* Extract the indices, assume there can be duplicate entries */
12254dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
122677c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1227e4d965acSSatish Balay 
1228e4d965acSSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1229e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
123006763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
12314dcbc457SBarry Smith     }
123206763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
123306763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1234e4d965acSSatish Balay 
123504a348a9SBarry Smith     k = 0;
123604a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
123704a348a9SBarry Smith       n = isz;
123806763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1239e4d965acSSatish Balay         row   = nidx[k];
1240e4d965acSSatish Balay         start = ai[row];
1241e4d965acSSatish Balay         end   = ai[row+1];
124204a348a9SBarry Smith         for ( l = start; l<end ; l++){
12438a047759SSatish Balay           val = aj[l] + shift;
124406763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1245e4d965acSSatish Balay         }
1246e4d965acSSatish Balay       }
1247e4d965acSSatish Balay     }
1248e4d965acSSatish Balay     ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1249e4d965acSSatish Balay   }
125004a348a9SBarry Smith   PetscFree(table);
125106763907SSatish Balay   PetscFree(nidx);
1252e4d965acSSatish Balay   return 0;
12534dcbc457SBarry Smith }
125417ab2063SBarry Smith 
1255682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1256682d7d0cSBarry Smith {
1257682d7d0cSBarry Smith   static int called = 0;
1258682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1259682d7d0cSBarry Smith 
1260682d7d0cSBarry Smith   if (called) return 0; else called = 1;
126177c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
126277c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
126377c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
126477c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
126577c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1266682d7d0cSBarry Smith #if defined(HAVE_ESSL)
126777c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1268682d7d0cSBarry Smith #endif
1269682d7d0cSBarry Smith   return 0;
1270682d7d0cSBarry Smith }
1271205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1272682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
127317ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
127417ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1275416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1276416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
127717ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
127817ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
127917ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
128017ab2063SBarry Smith        MatRelax_SeqAIJ,
128117ab2063SBarry Smith        MatTranspose_SeqAIJ,
12827264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1283f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
128417ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
128517ab2063SBarry Smith        MatCompress_SeqAIJ,
128617ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
128717ab2063SBarry Smith        MatGetReordering_SeqAIJ,
128817ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
128917ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
129017ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
129117ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1292416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
12933d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1294cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
12957eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1296682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1297f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
1298f0b747eeSBarry Smith        MatScale_SeqAIJ};
129917ab2063SBarry Smith 
130017ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
130117ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
130217ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
130317ab2063SBarry Smith 
130417ab2063SBarry Smith /*@C
1305682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
13060d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
13076e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
13080d15e28bSLois Curfman McInnes    (or nzz).  By setting these parameters accurately, performance can be
13090d15e28bSLois Curfman McInnes    increased by more than a factor of 50.
131017ab2063SBarry Smith 
131117ab2063SBarry Smith    Input Parameters:
131217ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
131317ab2063SBarry Smith .  m - number of rows
131417ab2063SBarry Smith .  n - number of columns
131517ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
131617ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
131717ab2063SBarry Smith 
131817ab2063SBarry Smith    Output Parameter:
1319416022c9SBarry Smith .  A - the matrix
132017ab2063SBarry Smith 
132117ab2063SBarry Smith    Notes:
132217ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
132317ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
13240002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
13250002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
132617ab2063SBarry Smith 
132717ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1328a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
13290d15e28bSLois Curfman McInnes    allocation.  For additional details, see the users manual chapter on
13300d15e28bSLois Curfman McInnes    matrices and the file $(PETSC_DIR)/Performance.
133117ab2063SBarry Smith 
1332682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1333682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1334682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
13356c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
13366c7ebb05SLois Curfman McInnes 
13376c7ebb05SLois Curfman McInnes    Options Database Keys:
13386c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
13396e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
13406e62573dSLois Curfman McInnes $        (max limit=5)
13416e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
13426e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
13436e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
134417ab2063SBarry Smith 
134517ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
134617ab2063SBarry Smith @*/
1347416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
134817ab2063SBarry Smith {
1349416022c9SBarry Smith   Mat        B;
1350416022c9SBarry Smith   Mat_SeqAIJ *b;
135169957df2SSatish Balay   int        i,len,ierr, flg;
1352d5d45c9bSBarry Smith 
1353416022c9SBarry Smith   *A      = 0;
13540452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1355416022c9SBarry Smith   PLogObjectCreate(B);
13560452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1357416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1358416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1359416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1360416022c9SBarry Smith   B->factor           = 0;
1361416022c9SBarry Smith   B->lupivotthreshold = 1.0;
13627a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
136369957df2SSatish Balay                           &flg); CHKERRQ(ierr);
13647a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
13657a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
13667a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1367416022c9SBarry Smith   b->row              = 0;
1368416022c9SBarry Smith   b->col              = 0;
1369416022c9SBarry Smith   b->indexshift       = 0;
1370b810aeb4SBarry Smith   b->reallocs         = 0;
137169957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
137269957df2SSatish Balay   if (flg) b->indexshift = -1;
137317ab2063SBarry Smith 
1374416022c9SBarry Smith   b->m       = m;
1375416022c9SBarry Smith   b->n       = n;
13760452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1377b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
13787b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
13797b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1380416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
138117ab2063SBarry Smith     nz = nz*m;
138217ab2063SBarry Smith   }
138317ab2063SBarry Smith   else {
138417ab2063SBarry Smith     nz = 0;
1385416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
138617ab2063SBarry Smith   }
138717ab2063SBarry Smith 
138817ab2063SBarry Smith   /* allocate the matrix space */
1389416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
13900452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1391416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1392cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1393416022c9SBarry Smith   b->i  = b->j + nz;
1394416022c9SBarry Smith   b->singlemalloc = 1;
139517ab2063SBarry Smith 
1396416022c9SBarry Smith   b->i[0] = -b->indexshift;
139717ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1398416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
139917ab2063SBarry Smith   }
140017ab2063SBarry Smith 
1401416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
14020452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1403416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1404416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
140517ab2063SBarry Smith 
1406416022c9SBarry Smith   b->nz               = 0;
1407416022c9SBarry Smith   b->maxnz            = nz;
1408416022c9SBarry Smith   b->sorted           = 0;
1409416022c9SBarry Smith   b->roworiented      = 1;
1410416022c9SBarry Smith   b->nonew            = 0;
1411416022c9SBarry Smith   b->diag             = 0;
1412416022c9SBarry Smith   b->solve_work       = 0;
141371bd300dSLois Curfman McInnes   b->spptr            = 0;
1414754ec7b1SSatish Balay   b->inode.node_count = 0;
1415754ec7b1SSatish Balay   b->inode.size       = 0;
14166c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
14176c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
141817ab2063SBarry Smith 
1419416022c9SBarry Smith   *A = B;
14204b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
14214b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
142269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
142369957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
14244b14c69eSBarry Smith #endif
142569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
142669957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
142769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
142869957df2SSatish Balay   if (flg) {
1429416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1430416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
143117ab2063SBarry Smith   }
143269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
143369957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
143417ab2063SBarry Smith   return 0;
143517ab2063SBarry Smith }
143617ab2063SBarry Smith 
14373d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
143817ab2063SBarry Smith {
1439416022c9SBarry Smith   Mat        C;
1440416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
144108480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
144217ab2063SBarry Smith 
14434043dd9cSLois Curfman McInnes   *B = 0;
14440452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1445416022c9SBarry Smith   PLogObjectCreate(C);
14460452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
144741c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1448416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1449416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1450416022c9SBarry Smith   C->factor     = A->factor;
1451416022c9SBarry Smith   c->row        = 0;
1452416022c9SBarry Smith   c->col        = 0;
1453416022c9SBarry Smith   c->indexshift = shift;
1454c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
145517ab2063SBarry Smith 
1456416022c9SBarry Smith   c->m          = a->m;
1457416022c9SBarry Smith   c->n          = a->n;
145817ab2063SBarry Smith 
14590452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
14600452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
146117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1462416022c9SBarry Smith     c->imax[i] = a->imax[i];
1463416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
146417ab2063SBarry Smith   }
146517ab2063SBarry Smith 
146617ab2063SBarry Smith   /* allocate the matrix space */
1467416022c9SBarry Smith   c->singlemalloc = 1;
1468416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
14690452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1470416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1471416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1472416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
147317ab2063SBarry Smith   if (m > 0) {
1474416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
147508480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1476416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
147717ab2063SBarry Smith     }
147808480c60SBarry Smith   }
147917ab2063SBarry Smith 
1480416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1481416022c9SBarry Smith   c->sorted      = a->sorted;
1482416022c9SBarry Smith   c->roworiented = a->roworiented;
1483416022c9SBarry Smith   c->nonew       = a->nonew;
14847a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
148517ab2063SBarry Smith 
1486416022c9SBarry Smith   if (a->diag) {
14870452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1488416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
148917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1490416022c9SBarry Smith       c->diag[i] = a->diag[i];
149117ab2063SBarry Smith     }
149217ab2063SBarry Smith   }
1493416022c9SBarry Smith   else c->diag          = 0;
14946c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
14956c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1496754ec7b1SSatish Balay   if (a->inode.size){
1497754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1498754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1499754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1500754ec7b1SSatish Balay   } else {
1501754ec7b1SSatish Balay     c->inode.size       = 0;
1502754ec7b1SSatish Balay     c->inode.node_count = 0;
1503754ec7b1SSatish Balay   }
1504416022c9SBarry Smith   c->nz                 = a->nz;
1505416022c9SBarry Smith   c->maxnz              = a->maxnz;
1506416022c9SBarry Smith   c->solve_work         = 0;
150776dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1508754ec7b1SSatish Balay 
1509416022c9SBarry Smith   *B = C;
151017ab2063SBarry Smith   return 0;
151117ab2063SBarry Smith }
151217ab2063SBarry Smith 
151319bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
151417ab2063SBarry Smith {
1515416022c9SBarry Smith   Mat_SeqAIJ   *a;
1516416022c9SBarry Smith   Mat          B;
151717699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1518bcd2baecSBarry Smith   MPI_Comm     comm;
151917ab2063SBarry Smith 
152019bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
152117699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
152217699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
152390ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
152477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
152548d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
152617ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
152717ab2063SBarry Smith 
152817ab2063SBarry Smith   /* read in row lengths */
15290452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
153077c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
153117ab2063SBarry Smith 
153217ab2063SBarry Smith   /* create our matrix */
1533416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1534416022c9SBarry Smith   B = *A;
1535416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1536416022c9SBarry Smith   shift = a->indexshift;
153717ab2063SBarry Smith 
153817ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
153977c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
154017ab2063SBarry Smith   if (shift) {
154117ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1542416022c9SBarry Smith       a->j[i] += 1;
154317ab2063SBarry Smith     }
154417ab2063SBarry Smith   }
154517ab2063SBarry Smith 
154617ab2063SBarry Smith   /* read in nonzero values */
154777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
154817ab2063SBarry Smith 
154917ab2063SBarry Smith   /* set matrix "i" values */
1550416022c9SBarry Smith   a->i[0] = -shift;
155117ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1552416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1553416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
155417ab2063SBarry Smith   }
15550452661fSBarry Smith   PetscFree(rowlengths);
155617ab2063SBarry Smith 
1557416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1558416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
155917ab2063SBarry Smith   return 0;
156017ab2063SBarry Smith }
156117ab2063SBarry Smith 
1562686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
15637264ac53SSatish Balay {
15647264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
15657264ac53SSatish Balay 
1566bcd2baecSBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type");
15677264ac53SSatish Balay 
15687264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
15697264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1570bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
157177c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1572bcd2baecSBarry Smith   }
15737264ac53SSatish Balay 
15747264ac53SSatish Balay   /* if the a->i are the same */
1575bcd2baecSBarry Smith   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
157677c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
15777264ac53SSatish Balay   }
15787264ac53SSatish Balay 
15797264ac53SSatish Balay   /* if a->j are the same */
1580bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
158177c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1582bcd2baecSBarry Smith   }
1583bcd2baecSBarry Smith 
1584bcd2baecSBarry Smith   /* if a->a are the same */
158519bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
158677c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
15877264ac53SSatish Balay   }
158877c4ece6SBarry Smith   *flg = PETSC_TRUE;
15897264ac53SSatish Balay   return 0;
15907264ac53SSatish Balay 
15917264ac53SSatish Balay }
1592