xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 95e01e2f099a96a22068d8621ba47172d7865843)
117ab2063SBarry Smith #ifndef lint
2*95e01e2fSLois Curfman McInnes static char vcid[] = "$Id: aij.c,v 1.150 1996/02/16 20:05:38 balay Exp curfman $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5d5d45c9bSBarry Smith /*
6d5d45c9bSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7d5d45c9bSBarry Smith   matrix storage format.
8d5d45c9bSBarry Smith */
917ab2063SBarry Smith #include "aij.h"
1017ab2063SBarry Smith #include "vec/vecimpl.h"
1117ab2063SBarry Smith #include "inline/spops.h"
12e4d965acSSatish Balay #include "petsc.h"
1306763907SSatish Balay #include "inline/bitarray.h"
1417ab2063SBarry Smith 
157a743949SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int**,int**);
1617ab2063SBarry Smith 
17416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
1817ab2063SBarry Smith {
19416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
20a2744918SBarry Smith   int        ierr, *ia, *ja,n,*idx,i;
213d54f168SSatish Balay   /*Viewer     V1, V2;*/
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 
427a743949SBarry Smith   ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,a->indexshift, &ia, &ja);CHKERRQ(ierr);
43416022c9SBarry Smith   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
44f31bba2fSSatish Balay 
450452661fSBarry Smith   PetscFree(ia); PetscFree(ja);
4617ab2063SBarry Smith   return 0;
4717ab2063SBarry Smith }
4817ab2063SBarry Smith 
49227d817aSBarry Smith #define CHUNKSIZE   15
5017ab2063SBarry Smith 
5117ab2063SBarry Smith /* This version has row oriented v  */
52416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
5317ab2063SBarry Smith {
54416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
55416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
564b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
57d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
58416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
5917ab2063SBarry Smith 
6017ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
61416022c9SBarry Smith     row  = im[k];
6217ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
63416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
6417ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
6517ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
66416022c9SBarry Smith     low = 0;
6717ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
68416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
69416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
704b0e389bSBarry Smith       col = in[l] - shift;
714b0e389bSBarry Smith       if (roworiented) {
724b0e389bSBarry Smith         value = *v++;
734b0e389bSBarry Smith       }
744b0e389bSBarry Smith       else {
754b0e389bSBarry Smith         value = v[k + l*m];
764b0e389bSBarry Smith       }
77416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
78416022c9SBarry Smith       while (high-low > 5) {
79416022c9SBarry Smith         t = (low+high)/2;
80416022c9SBarry Smith         if (rp[t] > col) high = t;
81416022c9SBarry Smith         else             low  = t;
8217ab2063SBarry Smith       }
83416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
8417ab2063SBarry Smith         if (rp[i] > col) break;
8517ab2063SBarry Smith         if (rp[i] == col) {
86416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
8717ab2063SBarry Smith           else                  ap[i] = value;
8817ab2063SBarry Smith           goto noinsert;
8917ab2063SBarry Smith         }
9017ab2063SBarry Smith       }
9117ab2063SBarry Smith       if (nonew) goto noinsert;
9217ab2063SBarry Smith       if (nrow >= rmax) {
9317ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
94416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
9517ab2063SBarry Smith         Scalar *new_a;
9617ab2063SBarry Smith 
9717ab2063SBarry Smith         /* malloc new storage space */
98416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
990452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
10017ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
10117ab2063SBarry Smith         new_i   = new_j + new_nz;
10217ab2063SBarry Smith 
10317ab2063SBarry Smith         /* copy over old data into new slots */
10417ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
105416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
106416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
107416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
108416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
10917ab2063SBarry Smith                                                            len*sizeof(int));
110416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
111416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
11217ab2063SBarry Smith                                                            len*sizeof(Scalar));
11317ab2063SBarry Smith         /* free up old matrix storage */
1140452661fSBarry Smith         PetscFree(a->a);
1150452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
116416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
117416022c9SBarry Smith         a->singlemalloc = 1;
11817ab2063SBarry Smith 
11917ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
120416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
121416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
122416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
123b810aeb4SBarry Smith         a->reallocs++;
12417ab2063SBarry Smith       }
125416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
126416022c9SBarry Smith       /* shift up all the later entries in this row */
127416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
12817ab2063SBarry Smith         rp[ii+1] = rp[ii];
12917ab2063SBarry Smith         ap[ii+1] = ap[ii];
13017ab2063SBarry Smith       }
13117ab2063SBarry Smith       rp[i] = col;
13217ab2063SBarry Smith       ap[i] = value;
13317ab2063SBarry Smith       noinsert:;
134416022c9SBarry Smith       low = i + 1;
13517ab2063SBarry Smith     }
13617ab2063SBarry Smith     ailen[row] = nrow;
13717ab2063SBarry Smith   }
13817ab2063SBarry Smith   return 0;
13917ab2063SBarry Smith }
14017ab2063SBarry Smith 
1417eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1427eb43aa7SLois Curfman McInnes {
1437eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
144b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1457eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1467eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
1477eb43aa7SLois Curfman McInnes 
1487eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
1497eb43aa7SLois Curfman McInnes     row  = im[k];
1507eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
1517eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
1527eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
1537eb43aa7SLois Curfman McInnes     nrow = ailen[row];
1547eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
1557eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
1567eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
1577eb43aa7SLois Curfman McInnes       col = in[l] - shift;
1587eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
1597eb43aa7SLois Curfman McInnes       while (high-low > 5) {
1607eb43aa7SLois Curfman McInnes         t = (low+high)/2;
1617eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
1627eb43aa7SLois Curfman McInnes         else             low  = t;
1637eb43aa7SLois Curfman McInnes       }
1647eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
1657eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
1667eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
167b49de8d1SLois Curfman McInnes           *v++ = ap[i];
1687eb43aa7SLois Curfman McInnes           goto finished;
1697eb43aa7SLois Curfman McInnes         }
1707eb43aa7SLois Curfman McInnes       }
171b49de8d1SLois Curfman McInnes       *v++ = zero;
1727eb43aa7SLois Curfman McInnes       finished:;
1737eb43aa7SLois Curfman McInnes     }
1747eb43aa7SLois Curfman McInnes   }
1757eb43aa7SLois Curfman McInnes   return 0;
1767eb43aa7SLois Curfman McInnes }
1777eb43aa7SLois Curfman McInnes 
17817ab2063SBarry Smith #include "draw.h"
17917ab2063SBarry Smith #include "pinclude/pviewer.h"
180416022c9SBarry Smith #include "sysio.h"
18117ab2063SBarry Smith 
182416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
18317ab2063SBarry Smith {
184416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
185416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
18617ab2063SBarry Smith 
187416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1880452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
189416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
190416022c9SBarry Smith   col_lens[1] = a->m;
191416022c9SBarry Smith   col_lens[2] = a->n;
192416022c9SBarry Smith   col_lens[3] = a->nz;
193416022c9SBarry Smith 
194416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
195416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
196416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
19717ab2063SBarry Smith   }
198416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
1990452661fSBarry Smith   PetscFree(col_lens);
200416022c9SBarry Smith 
201416022c9SBarry Smith   /* store column indices (zero start index) */
202416022c9SBarry Smith   if (a->indexshift) {
203416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
20417ab2063SBarry Smith   }
205416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
206416022c9SBarry Smith   if (a->indexshift) {
207416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
20817ab2063SBarry Smith   }
209416022c9SBarry Smith 
210416022c9SBarry Smith   /* store nonzero values */
211416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
21217ab2063SBarry Smith   return 0;
21317ab2063SBarry Smith }
214416022c9SBarry Smith 
215416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
216416022c9SBarry Smith {
217416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
218*95e01e2fSLois Curfman McInnes   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg;
21917ab2063SBarry Smith   FILE        *fd;
22017ab2063SBarry Smith   char        *outputname;
22117ab2063SBarry Smith 
222227d817aSBarry Smith   ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr);
223416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
224416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
22517ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
226*95e01e2fSLois Curfman McInnes     return 0;
227*95e01e2fSLois Curfman McInnes   }
228*95e01e2fSLois Curfman McInnes   else if (format == FILE_FORMAT_INFO_DETAILED) {
229*95e01e2fSLois Curfman McInnes     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr);
230*95e01e2fSLois Curfman McInnes     if (flg) fprintf(fd,"  not using I-node routines\n");
231*95e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
232*95e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
23317ab2063SBarry Smith   }
23417ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
23517ab2063SBarry Smith     int nz, nzalloc, mem;
236416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
237416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
23817ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
23917ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
24017ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
24117ab2063SBarry Smith 
24217ab2063SBarry Smith     for (i=0; i<m; i++) {
243416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
24417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
2457a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]),
246416022c9SBarry Smith                    imag(a->a[j]));
24717ab2063SBarry Smith #else
2487a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
24917ab2063SBarry Smith #endif
25017ab2063SBarry Smith       }
25117ab2063SBarry Smith     }
25217ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
25317ab2063SBarry Smith   }
25417ab2063SBarry Smith   else {
25517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
25617ab2063SBarry Smith       fprintf(fd,"row %d:",i);
257416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
25817ab2063SBarry Smith #if defined(PETSC_COMPLEX)
259416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
260416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
26117ab2063SBarry Smith         }
26217ab2063SBarry Smith         else {
263416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
26417ab2063SBarry Smith         }
26517ab2063SBarry Smith #else
266416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
26717ab2063SBarry Smith #endif
26817ab2063SBarry Smith       }
26917ab2063SBarry Smith       fprintf(fd,"\n");
27017ab2063SBarry Smith     }
27117ab2063SBarry Smith   }
27217ab2063SBarry Smith   fflush(fd);
273416022c9SBarry Smith   return 0;
274416022c9SBarry Smith }
275416022c9SBarry Smith 
276416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
277416022c9SBarry Smith {
278416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
279cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
280cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
281d7e8b826SBarry Smith   Draw     draw = (Draw) viewer;
282cddf8d76SBarry Smith   DrawButton  button;
283cddf8d76SBarry Smith 
284416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
285416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
286416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
287416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
288cddf8d76SBarry Smith   color = DRAW_BLUE;
289416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
290cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
291416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
292cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
293cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
294cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
295cddf8d76SBarry Smith #else
296cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
297cddf8d76SBarry Smith #endif
298cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
299cddf8d76SBarry Smith     }
300cddf8d76SBarry Smith   }
301cddf8d76SBarry Smith   color = DRAW_CYAN;
302cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
303cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
304cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
305cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
306cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
307cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
308cddf8d76SBarry Smith     }
309cddf8d76SBarry Smith   }
310cddf8d76SBarry Smith   color = DRAW_RED;
311cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
312cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
313cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
314cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
315cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
316cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
317cddf8d76SBarry Smith #else
318cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
319cddf8d76SBarry Smith #endif
320cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
321416022c9SBarry Smith     }
322416022c9SBarry Smith   }
323416022c9SBarry Smith   DrawFlush(draw);
324cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
325cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
326cddf8d76SBarry Smith 
327cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
328cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
329cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
330cddf8d76SBarry Smith     DrawClear(draw);
331cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
332cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
333cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
334cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
335cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
336cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
337cddf8d76SBarry Smith     w *= scale; h *= scale;
338cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
339cddf8d76SBarry Smith     color = DRAW_BLUE;
340cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
341cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
342cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
343cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
344cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
345cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
346cddf8d76SBarry Smith #else
347cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
348cddf8d76SBarry Smith #endif
349cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
350cddf8d76SBarry Smith       }
351cddf8d76SBarry Smith     }
352cddf8d76SBarry Smith     color = DRAW_CYAN;
353cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
354cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
355cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
356cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
357cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
358cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
359cddf8d76SBarry Smith       }
360cddf8d76SBarry Smith     }
361cddf8d76SBarry Smith     color = DRAW_RED;
362cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
363cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
364cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
365cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
366cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
367cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
368cddf8d76SBarry Smith #else
369cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
370cddf8d76SBarry Smith #endif
371cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
372cddf8d76SBarry Smith       }
373cddf8d76SBarry Smith     }
374cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
375cddf8d76SBarry Smith   }
376416022c9SBarry Smith   return 0;
377416022c9SBarry Smith }
378416022c9SBarry Smith 
379416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
380416022c9SBarry Smith {
381416022c9SBarry Smith   Mat         A = (Mat) obj;
382416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
383416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
384416022c9SBarry Smith 
385416022c9SBarry Smith   if (!viewer) {
386416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
387416022c9SBarry Smith   }
388416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
389416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
390416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
391416022c9SBarry Smith     }
392416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
393416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
394416022c9SBarry Smith     }
395416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
396416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
397416022c9SBarry Smith     }
398416022c9SBarry Smith   }
399416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
400416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
401416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
40217ab2063SBarry Smith   }
40317ab2063SBarry Smith   return 0;
40417ab2063SBarry Smith }
405c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
406416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
40717ab2063SBarry Smith {
408416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
40941c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
410416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
411416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
41217ab2063SBarry Smith 
41317ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
41417ab2063SBarry Smith 
41517ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
416416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
41717ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
41817ab2063SBarry Smith     if (fshift) {
419416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
42017ab2063SBarry Smith       N = ailen[i];
42117ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
42217ab2063SBarry Smith         ip[j-fshift] = ip[j];
42317ab2063SBarry Smith         ap[j-fshift] = ap[j];
42417ab2063SBarry Smith       }
42517ab2063SBarry Smith     }
42617ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
42717ab2063SBarry Smith   }
42817ab2063SBarry Smith   if (m) {
42917ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
43017ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
43117ab2063SBarry Smith   }
43217ab2063SBarry Smith   /* reset ilen and imax for each row */
43317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
43417ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
43517ab2063SBarry Smith   }
436416022c9SBarry Smith   a->nz = ai[m] + shift;
43717ab2063SBarry Smith 
43817ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
439416022c9SBarry Smith   if (fshift && a->diag) {
4400452661fSBarry Smith     PetscFree(a->diag);
441416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
442416022c9SBarry Smith     a->diag = 0;
44317ab2063SBarry Smith   }
444b810aeb4SBarry Smith   PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n",
445b810aeb4SBarry Smith            fshift,a->nz,m);
446b810aeb4SBarry Smith   PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n",
447b810aeb4SBarry Smith            a->reallocs);
44876dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
44941c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
45017ab2063SBarry Smith   return 0;
45117ab2063SBarry Smith }
45217ab2063SBarry Smith 
453416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
45417ab2063SBarry Smith {
455416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
456cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
45717ab2063SBarry Smith   return 0;
45817ab2063SBarry Smith }
459416022c9SBarry Smith 
46017ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
46117ab2063SBarry Smith {
462416022c9SBarry Smith   Mat        A  = (Mat) obj;
463416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
464d5d45c9bSBarry Smith 
46517ab2063SBarry Smith #if defined(PETSC_LOG)
466416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
46717ab2063SBarry Smith #endif
4680452661fSBarry Smith   PetscFree(a->a);
4690452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4700452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
4710452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4720452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
4730452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
47476dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
4750452661fSBarry Smith   PetscFree(a);
476416022c9SBarry Smith   PLogObjectDestroy(A);
4770452661fSBarry Smith   PetscHeaderDestroy(A);
47817ab2063SBarry Smith   return 0;
47917ab2063SBarry Smith }
48017ab2063SBarry Smith 
481416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
48217ab2063SBarry Smith {
48317ab2063SBarry Smith   return 0;
48417ab2063SBarry Smith }
48517ab2063SBarry Smith 
486416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
48717ab2063SBarry Smith {
488416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
489416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
4904b0e389bSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
491416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
492416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
493416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
494e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
495e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
496e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
497e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
498df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
499df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
5004d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
5016c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
5026c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
5036c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
5046c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
5056c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
506e2f28af5SBarry Smith   else
5074d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
50817ab2063SBarry Smith   return 0;
50917ab2063SBarry Smith }
51017ab2063SBarry Smith 
511416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
51217ab2063SBarry Smith {
513416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
514416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
51517ab2063SBarry Smith   Scalar     *x, zero = 0.0;
51617ab2063SBarry Smith 
51717ab2063SBarry Smith   VecSet(&zero,v);
51817ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
519416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
520416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
521416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
522416022c9SBarry Smith       if (a->j[j]+shift == i) {
523416022c9SBarry Smith         x[i] = a->a[j];
52417ab2063SBarry Smith         break;
52517ab2063SBarry Smith       }
52617ab2063SBarry Smith     }
52717ab2063SBarry Smith   }
52817ab2063SBarry Smith   return 0;
52917ab2063SBarry Smith }
53017ab2063SBarry Smith 
53117ab2063SBarry Smith /* -------------------------------------------------------*/
53217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
53317ab2063SBarry Smith /* -------------------------------------------------------*/
534416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
53517ab2063SBarry Smith {
536416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
53717ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
538416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
53917ab2063SBarry Smith 
54017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
541cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
54217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
54317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
544416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
545416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
546416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
54717ab2063SBarry Smith     alpha = x[i];
54817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
54917ab2063SBarry Smith   }
550416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
55117ab2063SBarry Smith   return 0;
55217ab2063SBarry Smith }
553d5d45c9bSBarry Smith 
554416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
55517ab2063SBarry Smith {
556416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
55717ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
558416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
55917ab2063SBarry Smith 
56017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
56117ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
56217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
56317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
564416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
565416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
566416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
56717ab2063SBarry Smith     alpha = x[i];
56817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
56917ab2063SBarry Smith   }
57017ab2063SBarry Smith   return 0;
57117ab2063SBarry Smith }
57217ab2063SBarry Smith 
573416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
57417ab2063SBarry Smith {
575416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
57617ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
577416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
57817ab2063SBarry Smith 
57917ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
58017ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
581416022c9SBarry Smith   idx  = a->j;
582416022c9SBarry Smith   v    = a->a;
583416022c9SBarry Smith   ii   = a->i;
58417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
585416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
58617ab2063SBarry Smith     sum  = 0.0;
58717ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
58817ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
589416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
59017ab2063SBarry Smith     y[i] = sum;
59117ab2063SBarry Smith   }
592416022c9SBarry Smith   PLogFlops(2*a->nz - m);
59317ab2063SBarry Smith   return 0;
59417ab2063SBarry Smith }
59517ab2063SBarry Smith 
596416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
59717ab2063SBarry Smith {
598416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
59917ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
600cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
60117ab2063SBarry Smith 
60217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
60317ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
604cddf8d76SBarry Smith   idx  = a->j;
605cddf8d76SBarry Smith   v    = a->a;
606cddf8d76SBarry Smith   ii   = a->i;
60717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
608cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
60917ab2063SBarry Smith     sum  = y[i];
610cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
611cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
61217ab2063SBarry Smith     z[i] = sum;
61317ab2063SBarry Smith   }
614416022c9SBarry Smith   PLogFlops(2*a->nz);
61517ab2063SBarry Smith   return 0;
61617ab2063SBarry Smith }
61717ab2063SBarry Smith 
61817ab2063SBarry Smith /*
61917ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
62017ab2063SBarry Smith */
62117ab2063SBarry Smith 
622416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
62317ab2063SBarry Smith {
624416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
625416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
62617ab2063SBarry Smith 
6270452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
628416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
629416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
630416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
631416022c9SBarry Smith       if (a->j[j]+shift == i) {
63217ab2063SBarry Smith         diag[i] = j - shift;
63317ab2063SBarry Smith         break;
63417ab2063SBarry Smith       }
63517ab2063SBarry Smith     }
63617ab2063SBarry Smith   }
637416022c9SBarry Smith   a->diag = diag;
63817ab2063SBarry Smith   return 0;
63917ab2063SBarry Smith }
64017ab2063SBarry Smith 
641416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
64217ab2063SBarry Smith                            double fshift,int its,Vec xx)
64317ab2063SBarry Smith {
644416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
645416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
646d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
64717ab2063SBarry Smith 
64817ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
649416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
650416022c9SBarry Smith   diag = a->diag;
651416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
65217ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
65317ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
65417ab2063SBarry Smith     bs = b + shift;
65517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
656416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
657416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
658416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
659416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
66017ab2063SBarry Smith         sum  = b[i]*d/omega;
66117ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
66217ab2063SBarry Smith         x[i] = sum;
66317ab2063SBarry Smith     }
66417ab2063SBarry Smith     return 0;
66517ab2063SBarry Smith   }
66617ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
66717ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
66817ab2063SBarry Smith   }
669416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
67017ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
67117ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
67217ab2063SBarry Smith 
67317ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
67417ab2063SBarry Smith 
67517ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
67617ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
67717ab2063SBarry Smith     is the relaxation factor.
67817ab2063SBarry Smith     */
6790452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
68017ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
68117ab2063SBarry Smith 
68217ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
68317ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
684416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
685416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
686416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
687416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
68817ab2063SBarry Smith       sum  = b[i];
68917ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
69017ab2063SBarry Smith       x[i] = omega*(sum/d);
69117ab2063SBarry Smith     }
69217ab2063SBarry Smith 
69317ab2063SBarry Smith     /*  t = b - (2*E - D)x */
694416022c9SBarry Smith     v = a->a;
69517ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
69617ab2063SBarry Smith 
69717ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
698416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
699416022c9SBarry Smith     diag = a->diag;
70017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
701416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
702416022c9SBarry Smith       n    = diag[i] - a->i[i];
703416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
704416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
70517ab2063SBarry Smith       sum  = t[i];
70617ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
70717ab2063SBarry Smith       t[i] = omega*(sum/d);
70817ab2063SBarry Smith     }
70917ab2063SBarry Smith 
71017ab2063SBarry Smith     /*  x = x + t */
71117ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7120452661fSBarry Smith     PetscFree(t);
71317ab2063SBarry Smith     return 0;
71417ab2063SBarry Smith   }
71517ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
71617ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
71717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
718416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
719416022c9SBarry Smith         n    = diag[i] - a->i[i];
720416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
721416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
72217ab2063SBarry Smith         sum  = b[i];
72317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
72417ab2063SBarry Smith         x[i] = omega*(sum/d);
72517ab2063SBarry Smith       }
72617ab2063SBarry Smith       xb = x;
72717ab2063SBarry Smith     }
72817ab2063SBarry Smith     else xb = b;
72917ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
73017ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
73117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
732416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
73317ab2063SBarry Smith       }
73417ab2063SBarry Smith     }
73517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
73617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
737416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
738416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
739416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
740416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
74117ab2063SBarry Smith         sum  = xb[i];
74217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
74317ab2063SBarry Smith         x[i] = omega*(sum/d);
74417ab2063SBarry Smith       }
74517ab2063SBarry Smith     }
74617ab2063SBarry Smith     its--;
74717ab2063SBarry Smith   }
74817ab2063SBarry Smith   while (its--) {
74917ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
75017ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
751416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
752416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
753416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
754416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
75517ab2063SBarry Smith         sum  = b[i];
75617ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
75717ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
75817ab2063SBarry Smith       }
75917ab2063SBarry Smith     }
76017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
76117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
762416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
763416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
764416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
765416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
76617ab2063SBarry Smith         sum  = b[i];
76717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
76817ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
76917ab2063SBarry Smith       }
77017ab2063SBarry Smith     }
77117ab2063SBarry Smith   }
77217ab2063SBarry Smith   return 0;
77317ab2063SBarry Smith }
77417ab2063SBarry Smith 
775d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
77617ab2063SBarry Smith {
777416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
778416022c9SBarry Smith   *nz      = a->nz;
779416022c9SBarry Smith   *nzalloc = a->maxnz;
780416022c9SBarry Smith   *mem     = (int)A->mem;
78117ab2063SBarry Smith   return 0;
78217ab2063SBarry Smith }
78317ab2063SBarry Smith 
78417ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
78517ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
78617ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
78717ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
78817ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
78917ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
79017ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
79117ab2063SBarry Smith 
79217ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
79317ab2063SBarry Smith {
794416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
795416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
79617ab2063SBarry Smith 
79717ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
79817ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
79917ab2063SBarry Smith   if (diag) {
80017ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
801416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
802416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
803416022c9SBarry Smith         a->ilen[rows[i]] = 1;
804416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
805416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
80617ab2063SBarry Smith       }
80717ab2063SBarry Smith       else {
80817ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
80917ab2063SBarry Smith         CHKERRQ(ierr);
81017ab2063SBarry Smith       }
81117ab2063SBarry Smith     }
81217ab2063SBarry Smith   }
81317ab2063SBarry Smith   else {
81417ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
815416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
816416022c9SBarry Smith       a->ilen[rows[i]] = 0;
81717ab2063SBarry Smith     }
81817ab2063SBarry Smith   }
81917ab2063SBarry Smith   ISRestoreIndices(is,&rows);
82017ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82117ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
82217ab2063SBarry Smith   return 0;
82317ab2063SBarry Smith }
82417ab2063SBarry Smith 
825416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
82617ab2063SBarry Smith {
827416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
828416022c9SBarry Smith   *m = a->m; *n = a->n;
82917ab2063SBarry Smith   return 0;
83017ab2063SBarry Smith }
83117ab2063SBarry Smith 
832416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
83317ab2063SBarry Smith {
834416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
835416022c9SBarry Smith   *m = 0; *n = a->m;
83617ab2063SBarry Smith   return 0;
83717ab2063SBarry Smith }
838416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
83917ab2063SBarry Smith {
840416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
841c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
84217ab2063SBarry Smith 
843416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
84417ab2063SBarry Smith 
845416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
846416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
84717ab2063SBarry Smith   if (idx) {
84817ab2063SBarry Smith     if (*nz) {
849416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
8500452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
85117ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
85217ab2063SBarry Smith     }
85317ab2063SBarry Smith     else *idx = 0;
85417ab2063SBarry Smith   }
85517ab2063SBarry Smith   return 0;
85617ab2063SBarry Smith }
85717ab2063SBarry Smith 
858416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
85917ab2063SBarry Smith {
8600452661fSBarry Smith   if (idx) {if (*idx) PetscFree(*idx);}
86117ab2063SBarry Smith   return 0;
86217ab2063SBarry Smith }
86317ab2063SBarry Smith 
864cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
86517ab2063SBarry Smith {
866416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
867416022c9SBarry Smith   Scalar     *v = a->a;
86817ab2063SBarry Smith   double     sum = 0.0;
869416022c9SBarry Smith   int        i, j,shift = a->indexshift;
87017ab2063SBarry Smith 
87117ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
872416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
87317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
87417ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
87517ab2063SBarry Smith #else
87617ab2063SBarry Smith       sum += (*v)*(*v); v++;
87717ab2063SBarry Smith #endif
87817ab2063SBarry Smith     }
87917ab2063SBarry Smith     *norm = sqrt(sum);
88017ab2063SBarry Smith   }
88117ab2063SBarry Smith   else if (type == NORM_1) {
88217ab2063SBarry Smith     double *tmp;
883416022c9SBarry Smith     int    *jj = a->j;
8840452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
885cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
88617ab2063SBarry Smith     *norm = 0.0;
887416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
888a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
88917ab2063SBarry Smith     }
890416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
89117ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
89217ab2063SBarry Smith     }
8930452661fSBarry Smith     PetscFree(tmp);
89417ab2063SBarry Smith   }
89517ab2063SBarry Smith   else if (type == NORM_INFINITY) {
89617ab2063SBarry Smith     *norm = 0.0;
897416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
898416022c9SBarry Smith       v = a->a + a->i[j] + shift;
89917ab2063SBarry Smith       sum = 0.0;
900416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
901cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
90217ab2063SBarry Smith       }
90317ab2063SBarry Smith       if (sum > *norm) *norm = sum;
90417ab2063SBarry Smith     }
90517ab2063SBarry Smith   }
90617ab2063SBarry Smith   else {
90748d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
90817ab2063SBarry Smith   }
90917ab2063SBarry Smith   return 0;
91017ab2063SBarry Smith }
91117ab2063SBarry Smith 
912416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
91317ab2063SBarry Smith {
914416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
915416022c9SBarry Smith   Mat        C;
916416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
917416022c9SBarry Smith   int        shift = a->indexshift;
918d5d45c9bSBarry Smith   Scalar     *array = a->a;
91917ab2063SBarry Smith 
9203638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
9213638b69dSLois Curfman McInnes     SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place");
9220452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
923cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
92417ab2063SBarry Smith   if (shift) {
92517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
92617ab2063SBarry Smith   }
92717ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
928416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
9290452661fSBarry Smith   PetscFree(col);
93017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
93117ab2063SBarry Smith     len = ai[i+1]-ai[i];
932416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
93317ab2063SBarry Smith     array += len; aj += len;
93417ab2063SBarry Smith   }
93517ab2063SBarry Smith   if (shift) {
93617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
93717ab2063SBarry Smith   }
93817ab2063SBarry Smith 
939416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
940416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
94117ab2063SBarry Smith 
9423638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
943416022c9SBarry Smith     *B = C;
94417ab2063SBarry Smith   } else {
945416022c9SBarry Smith     /* This isn't really an in-place transpose */
9460452661fSBarry Smith     PetscFree(a->a);
9470452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
9480452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
9490452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
9500452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
9510452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
9520452661fSBarry Smith     PetscFree(a);
953416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
9540452661fSBarry Smith     PetscHeaderDestroy(C);
95517ab2063SBarry Smith   }
95617ab2063SBarry Smith   return 0;
95717ab2063SBarry Smith }
95817ab2063SBarry Smith 
959f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
96017ab2063SBarry Smith {
961416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
96217ab2063SBarry Smith   Scalar     *l,*r,x,*v;
963d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
96417ab2063SBarry Smith 
96517ab2063SBarry Smith   if (ll) {
96617ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
967f0b747eeSBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
968416022c9SBarry Smith     v = a->a;
96917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
97017ab2063SBarry Smith       x = l[i];
971416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
97217ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
97317ab2063SBarry Smith     }
97417ab2063SBarry Smith   }
97517ab2063SBarry Smith   if (rr) {
97617ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
977f0b747eeSBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
978416022c9SBarry Smith     v = a->a; jj = a->j;
97917ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
98017ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
98117ab2063SBarry Smith     }
98217ab2063SBarry Smith   }
98317ab2063SBarry Smith   return 0;
98417ab2063SBarry Smith }
98517ab2063SBarry Smith 
986cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
98717ab2063SBarry Smith {
988db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
98902834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
990a2744918SBarry Smith   register int sum,lensi;
99102834360SBarry Smith   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
992a2744918SBarry Smith   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
993db02288aSLois Curfman McInnes   Scalar       *vwork,*a_new;
994416022c9SBarry Smith   Mat          C;
99517ab2063SBarry Smith 
99617ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
99717ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
99817ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
99917ab2063SBarry Smith 
10007264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
100102834360SBarry Smith     /* special case of contiguous rows */
10020452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
100302834360SBarry Smith     starts = lens + ncols;
100402834360SBarry Smith     /* loop over new rows determining lens and starting points */
100502834360SBarry Smith     for (i=0; i<nrows; i++) {
1006a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1007a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
100802834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1009d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
101002834360SBarry Smith           starts[i] = k;
101102834360SBarry Smith           break;
101202834360SBarry Smith 	}
101302834360SBarry Smith       }
1014a2744918SBarry Smith       sum = 0;
101502834360SBarry Smith       while (k < kend) {
1016d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1017a2744918SBarry Smith         sum++;
101802834360SBarry Smith       }
1019a2744918SBarry Smith       lens[i] = sum;
102002834360SBarry Smith     }
102102834360SBarry Smith     /* create submatrix */
1022cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
102308480c60SBarry Smith       int n_cols,n_rows;
102408480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
102508480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1026d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
102708480c60SBarry Smith       C = *B;
102808480c60SBarry Smith     }
102908480c60SBarry Smith     else {
103002834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
103108480c60SBarry Smith     }
1032db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1033db02288aSLois Curfman McInnes 
103402834360SBarry Smith     /* loop over rows inserting into submatrix */
1035db02288aSLois Curfman McInnes     a_new    = c->a;
1036db02288aSLois Curfman McInnes     j_new    = c->j;
1037db02288aSLois Curfman McInnes     i_new    = c->i;
1038db02288aSLois Curfman McInnes     i_new[0] = -shift;
103902834360SBarry Smith     for (i=0; i<nrows; i++) {
1040a2744918SBarry Smith       ii    = starts[i];
1041a2744918SBarry Smith       lensi = lens[i];
1042a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1043a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
104402834360SBarry Smith       }
1045a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1046a2744918SBarry Smith       a_new      += lensi;
1047a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1048a2744918SBarry Smith       c->ilen[i]  = lensi;
104902834360SBarry Smith     }
10500452661fSBarry Smith     PetscFree(lens);
105102834360SBarry Smith   }
105202834360SBarry Smith   else {
105302834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10547264ac53SSatish Balay     ierr = SYIsort(ncols,icol); CHKERRQ(ierr);
10550452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
105602834360SBarry Smith     ssmap = smap + shift;
10577eb43aa7SLois Curfman McInnes     cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork);
105802834360SBarry Smith     lens  = cwork + ncols;
10590452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1060cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
106117ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
106202834360SBarry Smith     /* determine lens of each row */
106302834360SBarry Smith     for (i=0; i<nrows; i++) {
1064d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
106502834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
106602834360SBarry Smith       lens[i] = 0;
106702834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1068d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
106902834360SBarry Smith           lens[i]++;
107002834360SBarry Smith         }
107102834360SBarry Smith       }
107202834360SBarry Smith     }
107317ab2063SBarry Smith     /* Create and fill new matrix */
1074a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
107508480c60SBarry Smith       int n_cols,n_rows;
107608480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
107708480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1078d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
107908480c60SBarry Smith       C = *B;
108008480c60SBarry Smith     }
108108480c60SBarry Smith     else {
108202834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
108308480c60SBarry Smith     }
108417ab2063SBarry Smith     for (i=0; i<nrows; i++) {
108517ab2063SBarry Smith       nznew  = 0;
1086d8ced48eSBarry Smith       kstart = ai[irow[i]]+shift;
1087416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
108817ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
108902834360SBarry Smith         if (ssmap[a->j[k]]) {
109002834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1091416022c9SBarry Smith           vwork[nznew++] = a->a[k];
109217ab2063SBarry Smith         }
109317ab2063SBarry Smith       }
1094416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
109517ab2063SBarry Smith     }
109602834360SBarry Smith     /* Free work space */
109702834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10980452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
109902834360SBarry Smith   }
1100416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1101416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
110217ab2063SBarry Smith 
110317ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1104416022c9SBarry Smith   *B = C;
110517ab2063SBarry Smith   return 0;
110617ab2063SBarry Smith }
110717ab2063SBarry Smith 
1108a871dcd8SBarry Smith /*
110963b91edcSBarry Smith      note: This can only work for identity for row and col. It would
111063b91edcSBarry Smith    be good to check this and otherwise generate an error.
1111a871dcd8SBarry Smith */
111263b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1113a871dcd8SBarry Smith {
111463b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
111508480c60SBarry Smith   int        ierr;
111663b91edcSBarry Smith   Mat        outA;
111763b91edcSBarry Smith 
1118a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1119a871dcd8SBarry Smith 
112063b91edcSBarry Smith   outA          = inA;
112163b91edcSBarry Smith   inA->factor   = FACTOR_LU;
112263b91edcSBarry Smith   a->row        = row;
112363b91edcSBarry Smith   a->col        = col;
112463b91edcSBarry Smith 
11250452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
112663b91edcSBarry Smith 
112708480c60SBarry Smith   if (!a->diag) {
112808480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
112963b91edcSBarry Smith   }
113063b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1131a871dcd8SBarry Smith   return 0;
1132a871dcd8SBarry Smith }
1133a871dcd8SBarry Smith 
1134f0b747eeSBarry Smith #include "pinclude/plapack.h"
1135f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1136f0b747eeSBarry Smith {
1137f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1138f0b747eeSBarry Smith   int        one = 1;
1139f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1140f0b747eeSBarry Smith   PLogFlops(a->nz);
1141f0b747eeSBarry Smith   return 0;
1142f0b747eeSBarry Smith }
1143f0b747eeSBarry Smith 
1144cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1145cddf8d76SBarry Smith                                     Mat **B)
1146cddf8d76SBarry Smith {
1147cddf8d76SBarry Smith   int ierr,i;
1148cddf8d76SBarry Smith 
1149cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
11500452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1151cddf8d76SBarry Smith   }
1152cddf8d76SBarry Smith 
1153cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1154cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1155cddf8d76SBarry Smith   }
1156cddf8d76SBarry Smith   return 0;
1157cddf8d76SBarry Smith }
1158cddf8d76SBarry Smith 
1159e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
11604dcbc457SBarry Smith {
1161e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
116206763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
11638a047759SSatish Balay   int        start, end, *ai, *aj;
116406763907SSatish Balay   char       *table;
11658a047759SSatish Balay   shift = a->indexshift;
1166e4d965acSSatish Balay   m     = a->m;
1167e4d965acSSatish Balay   ai    = a->i;
11688a047759SSatish Balay   aj    = a->j+shift;
11698a047759SSatish Balay 
1170e4d965acSSatish Balay   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
117106763907SSatish Balay 
117206763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
117306763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
117406763907SSatish Balay 
1175e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1176e4d965acSSatish Balay     /* Initialise the two local arrays */
1177e4d965acSSatish Balay     isz  = 0;
117806763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1179e4d965acSSatish Balay 
1180e4d965acSSatish Balay                 /* Extract the indices, assume there can be duplicate entries */
11814dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
1182e4d965acSSatish Balay     ierr = ISGetLocalSize(is[i],&n);  CHKERRQ(ierr);
1183e4d965acSSatish Balay 
1184e4d965acSSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1185e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
118606763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
11874dcbc457SBarry Smith     }
118806763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
118906763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1190e4d965acSSatish Balay 
119104a348a9SBarry Smith     k = 0;
119204a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
119304a348a9SBarry Smith       n = isz;
119406763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1195e4d965acSSatish Balay         row   = nidx[k];
1196e4d965acSSatish Balay         start = ai[row];
1197e4d965acSSatish Balay         end   = ai[row+1];
119804a348a9SBarry Smith         for ( l = start; l<end ; l++){
11998a047759SSatish Balay           val = aj[l] + shift;
120006763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1201e4d965acSSatish Balay         }
1202e4d965acSSatish Balay       }
1203e4d965acSSatish Balay     }
1204e4d965acSSatish Balay     ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1205e4d965acSSatish Balay   }
120604a348a9SBarry Smith   PetscFree(table);
120706763907SSatish Balay   PetscFree(nidx);
1208e4d965acSSatish Balay   return 0;
12094dcbc457SBarry Smith }
121017ab2063SBarry Smith 
1211682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1212682d7d0cSBarry Smith {
1213682d7d0cSBarry Smith   static int called = 0;
1214682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1215682d7d0cSBarry Smith 
1216682d7d0cSBarry Smith   if (called) return 0; else called = 1;
12172a7368beSLois Curfman McInnes   MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1218682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
12192a7368beSLois Curfman McInnes   MPIU_printf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
1220682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
1221682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1222682d7d0cSBarry Smith #if defined(HAVE_ESSL)
1223682d7d0cSBarry Smith   MPIU_printf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1224682d7d0cSBarry Smith #endif
1225682d7d0cSBarry Smith   return 0;
1226682d7d0cSBarry Smith }
12277264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg);
1228682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
122917ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
123017ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1231416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1232416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
123317ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
123417ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
123517ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
123617ab2063SBarry Smith        MatRelax_SeqAIJ,
123717ab2063SBarry Smith        MatTranspose_SeqAIJ,
12387264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1239f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
124017ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
124117ab2063SBarry Smith        MatCompress_SeqAIJ,
124217ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
124317ab2063SBarry Smith        MatGetReordering_SeqAIJ,
124417ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
124517ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
124617ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
124717ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1248416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
12493d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1250cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
12517eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1252682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1253f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
1254f0b747eeSBarry Smith        MatScale_SeqAIJ};
125517ab2063SBarry Smith 
125617ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
125717ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
125817ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
125917ab2063SBarry Smith 
126017ab2063SBarry Smith /*@C
1261682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
12620d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
12636e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
12640d15e28bSLois Curfman McInnes    (or nzz).  By setting these parameters accurately, performance can be
12650d15e28bSLois Curfman McInnes    increased by more than a factor of 50.
126617ab2063SBarry Smith 
126717ab2063SBarry Smith    Input Parameters:
126817ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
126917ab2063SBarry Smith .  m - number of rows
127017ab2063SBarry Smith .  n - number of columns
127117ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
127217ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
127317ab2063SBarry Smith 
127417ab2063SBarry Smith    Output Parameter:
1275416022c9SBarry Smith .  A - the matrix
127617ab2063SBarry Smith 
127717ab2063SBarry Smith    Notes:
127817ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
127917ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
12800002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
12810002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
128217ab2063SBarry Smith 
128317ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1284a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
12850d15e28bSLois Curfman McInnes    allocation.  For additional details, see the users manual chapter on
12860d15e28bSLois Curfman McInnes    matrices and the file $(PETSC_DIR)/Performance.
128717ab2063SBarry Smith 
1288682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1289682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1290682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
12916c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
12926c7ebb05SLois Curfman McInnes 
12936c7ebb05SLois Curfman McInnes    Options Database Keys:
12946c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
12956e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
12966e62573dSLois Curfman McInnes $        (max limit=5)
12976e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
12986e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
12996e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
130017ab2063SBarry Smith 
130117ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
130217ab2063SBarry Smith @*/
1303416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
130417ab2063SBarry Smith {
1305416022c9SBarry Smith   Mat        B;
1306416022c9SBarry Smith   Mat_SeqAIJ *b;
130769957df2SSatish Balay   int        i,len,ierr, flg;
1308d5d45c9bSBarry Smith 
1309416022c9SBarry Smith   *A      = 0;
13100452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1311416022c9SBarry Smith   PLogObjectCreate(B);
13120452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1313416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1314416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1315416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1316416022c9SBarry Smith   B->factor           = 0;
1317416022c9SBarry Smith   B->lupivotthreshold = 1.0;
13187a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
131969957df2SSatish Balay                           &flg); CHKERRQ(ierr);
13207a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
13217a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
13227a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1323416022c9SBarry Smith   b->row              = 0;
1324416022c9SBarry Smith   b->col              = 0;
1325416022c9SBarry Smith   b->indexshift       = 0;
1326b810aeb4SBarry Smith   b->reallocs         = 0;
132769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
132869957df2SSatish Balay   if (flg) b->indexshift = -1;
132917ab2063SBarry Smith 
1330416022c9SBarry Smith   b->m       = m;
1331416022c9SBarry Smith   b->n       = n;
13320452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1333b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
13347b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
13357b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1336416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
133717ab2063SBarry Smith     nz = nz*m;
133817ab2063SBarry Smith   }
133917ab2063SBarry Smith   else {
134017ab2063SBarry Smith     nz = 0;
1341416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
134217ab2063SBarry Smith   }
134317ab2063SBarry Smith 
134417ab2063SBarry Smith   /* allocate the matrix space */
1345416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
13460452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1347416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1348cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1349416022c9SBarry Smith   b->i  = b->j + nz;
1350416022c9SBarry Smith   b->singlemalloc = 1;
135117ab2063SBarry Smith 
1352416022c9SBarry Smith   b->i[0] = -b->indexshift;
135317ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1354416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
135517ab2063SBarry Smith   }
135617ab2063SBarry Smith 
1357416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
13580452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1359416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1360416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
136117ab2063SBarry Smith 
1362416022c9SBarry Smith   b->nz               = 0;
1363416022c9SBarry Smith   b->maxnz            = nz;
1364416022c9SBarry Smith   b->sorted           = 0;
1365416022c9SBarry Smith   b->roworiented      = 1;
1366416022c9SBarry Smith   b->nonew            = 0;
1367416022c9SBarry Smith   b->diag             = 0;
1368416022c9SBarry Smith   b->solve_work       = 0;
136971bd300dSLois Curfman McInnes   b->spptr            = 0;
1370754ec7b1SSatish Balay   b->inode.node_count = 0;
1371754ec7b1SSatish Balay   b->inode.size       = 0;
13726c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
13736c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
137417ab2063SBarry Smith 
1375416022c9SBarry Smith   *A = B;
13764b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
13774b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
137869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
137969957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
13804b14c69eSBarry Smith #endif
138169957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
138269957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
138369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
138469957df2SSatish Balay   if (flg) {
1385416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1386416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
138717ab2063SBarry Smith   }
138869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
138969957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
139017ab2063SBarry Smith   return 0;
139117ab2063SBarry Smith }
139217ab2063SBarry Smith 
13933d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
139417ab2063SBarry Smith {
1395416022c9SBarry Smith   Mat        C;
1396416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
139708480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
139817ab2063SBarry Smith 
13994043dd9cSLois Curfman McInnes   *B = 0;
14000452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1401416022c9SBarry Smith   PLogObjectCreate(C);
14020452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
140341c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1404416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1405416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1406416022c9SBarry Smith   C->factor     = A->factor;
1407416022c9SBarry Smith   c->row        = 0;
1408416022c9SBarry Smith   c->col        = 0;
1409416022c9SBarry Smith   c->indexshift = shift;
1410c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
141117ab2063SBarry Smith 
1412416022c9SBarry Smith   c->m          = a->m;
1413416022c9SBarry Smith   c->n          = a->n;
141417ab2063SBarry Smith 
14150452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
14160452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
141717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1418416022c9SBarry Smith     c->imax[i] = a->imax[i];
1419416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
142017ab2063SBarry Smith   }
142117ab2063SBarry Smith 
142217ab2063SBarry Smith   /* allocate the matrix space */
1423416022c9SBarry Smith   c->singlemalloc = 1;
1424416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
14250452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1426416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1427416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1428416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
142917ab2063SBarry Smith   if (m > 0) {
1430416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
143108480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1432416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
143317ab2063SBarry Smith     }
143408480c60SBarry Smith   }
143517ab2063SBarry Smith 
1436416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1437416022c9SBarry Smith   c->sorted      = a->sorted;
1438416022c9SBarry Smith   c->roworiented = a->roworiented;
1439416022c9SBarry Smith   c->nonew       = a->nonew;
14407a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
144117ab2063SBarry Smith 
1442416022c9SBarry Smith   if (a->diag) {
14430452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1444416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
144517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1446416022c9SBarry Smith       c->diag[i] = a->diag[i];
144717ab2063SBarry Smith     }
144817ab2063SBarry Smith   }
1449416022c9SBarry Smith   else c->diag          = 0;
14506c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
14516c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1452754ec7b1SSatish Balay   if (a->inode.size){
1453754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1454754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1455754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1456754ec7b1SSatish Balay   } else {
1457754ec7b1SSatish Balay     c->inode.size       = 0;
1458754ec7b1SSatish Balay     c->inode.node_count = 0;
1459754ec7b1SSatish Balay   }
1460416022c9SBarry Smith   c->nz                 = a->nz;
1461416022c9SBarry Smith   c->maxnz              = a->maxnz;
1462416022c9SBarry Smith   c->solve_work         = 0;
146376dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1464754ec7b1SSatish Balay 
1465416022c9SBarry Smith   *B = C;
146617ab2063SBarry Smith   return 0;
146717ab2063SBarry Smith }
146817ab2063SBarry Smith 
1469416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
147017ab2063SBarry Smith {
1471416022c9SBarry Smith   Mat_SeqAIJ   *a;
1472416022c9SBarry Smith   Mat          B;
147317699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
147417ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
147517ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
147617ab2063SBarry Smith 
147717699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
147817699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
147917ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1480416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
148148d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
148217ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
148317ab2063SBarry Smith 
148417ab2063SBarry Smith   /* read in row lengths */
14850452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1486416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
148717ab2063SBarry Smith 
148817ab2063SBarry Smith   /* create our matrix */
1489416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1490416022c9SBarry Smith   B = *A;
1491416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1492416022c9SBarry Smith   shift = a->indexshift;
149317ab2063SBarry Smith 
149417ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1495416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
149617ab2063SBarry Smith   if (shift) {
149717ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1498416022c9SBarry Smith       a->j[i] += 1;
149917ab2063SBarry Smith     }
150017ab2063SBarry Smith   }
150117ab2063SBarry Smith 
150217ab2063SBarry Smith   /* read in nonzero values */
1503416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
150417ab2063SBarry Smith 
150517ab2063SBarry Smith   /* set matrix "i" values */
1506416022c9SBarry Smith   a->i[0] = -shift;
150717ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1508416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1509416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
151017ab2063SBarry Smith   }
15110452661fSBarry Smith   PetscFree(rowlengths);
151217ab2063SBarry Smith 
1513416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1514416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
151517ab2063SBarry Smith   return 0;
151617ab2063SBarry Smith }
151717ab2063SBarry Smith 
151817ab2063SBarry Smith 
151917ab2063SBarry Smith 
15207264ac53SSatish Balay /*
15217264ac53SSatish Balay 
15227264ac53SSatish Balay   flg =0 if error;
15237264ac53SSatish Balay   flg =1 if both are equal;
15247264ac53SSatish Balay   */
15257264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg)
15267264ac53SSatish Balay {
15277264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
15287264ac53SSatish Balay 
15297264ac53SSatish Balay   if (B->type !=MATSEQAIJ)  SETERRQ(1,"MatEqual_SeqAIJ:Both matrices should be of type MATSEQAIJ");
15307264ac53SSatish Balay 
15317264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
15327264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
15337264ac53SSatish Balay       (a->indexshift != b->indexshift)) { *flg =0 ; return 0; }
15347264ac53SSatish Balay 
15357264ac53SSatish Balay   /* if the a->i are the same */
15367264ac53SSatish Balay   if(PetscMemcmp((char *)a->i, (char*)b->i, (a->n+1)*sizeof(int))) { *flg =0 ; return 0;}
15377264ac53SSatish Balay 
15387264ac53SSatish Balay   /* if a->j are the same */
15397264ac53SSatish Balay   if(PetscMemcmp((char *)a->j, (char*)b->j, (a->nz)*sizeof(int))) {
15407264ac53SSatish Balay     *flg =0 ; return 0;
15417264ac53SSatish Balay   }
15427264ac53SSatish Balay 
15437264ac53SSatish Balay   /* if a->j are the same */
15447264ac53SSatish Balay   if(PetscMemcmp((char *)a->a, (char*)b->a, (a->nz)*sizeof(Scalar))) {
15457264ac53SSatish Balay     *flg =0 ; return 0;
15467264ac53SSatish Balay   }
15477264ac53SSatish Balay   *flg =1 ;
15487264ac53SSatish Balay   return 0;
15497264ac53SSatish Balay 
15507264ac53SSatish Balay }
1551