xref: /petsc/src/mat/impls/aij/seq/aij.c (revision b49de8d1b1f74f27495c5d351d287e6bdc71c90c)
117ab2063SBarry Smith #ifndef lint
2*b49de8d1SLois Curfman McInnes static char vcid[] = "$Id: aij.c,v 1.121 1995/11/25 19:07:38 curfman 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"
1217ab2063SBarry Smith 
1317ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**);
1417ab2063SBarry Smith 
15416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
1617ab2063SBarry Smith {
17416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
18a2744918SBarry Smith   int        ierr, *ia, *ja,n,*idx,i;
193d54f168SSatish Balay   /*Viewer     V1, V2;*/
2017ab2063SBarry Smith 
21416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix");
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 
42416022c9SBarry Smith   ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr);
43416022c9SBarry Smith   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
4444f6d32bSSatish Balay /*  ISView(*rperm, STDOUT_VIEWER_SELF);*/
45f31bba2fSSatish Balay 
463d54f168SSatish Balay   /*ViewerFileOpenASCII(MPI_COMM_SELF,"row_is_orig", &V1);
47f31bba2fSSatish Balay   ViewerFileOpenASCII(MPI_COMM_SELF,"col_is_orig", &V2);
48f31bba2fSSatish Balay   ISView(*rperm,V1);
49f31bba2fSSatish Balay   ISView(*cperm,V2);
50f31bba2fSSatish Balay   ViewerDestroy(V1);
513d54f168SSatish Balay   ViewerDestroy(V2);*/
52f31bba2fSSatish Balay 
530452661fSBarry Smith   PetscFree(ia); PetscFree(ja);
5417ab2063SBarry Smith   return 0;
5517ab2063SBarry Smith }
5617ab2063SBarry Smith 
57416022c9SBarry Smith #define CHUNKSIZE   10
5817ab2063SBarry Smith 
5917ab2063SBarry Smith /* This version has row oriented v  */
60416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
6117ab2063SBarry Smith {
62416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
63416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
64416022c9SBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen;
65d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
66416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
6717ab2063SBarry Smith 
6817ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
69416022c9SBarry Smith     row  = im[k];
7017ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
71416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
7217ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
7317ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
74416022c9SBarry Smith     low = 0;
7517ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
76416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
77416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
78416022c9SBarry Smith       col = in[l] - shift; value = *v++;
79416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
80416022c9SBarry Smith       while (high-low > 5) {
81416022c9SBarry Smith         t = (low+high)/2;
82416022c9SBarry Smith         if (rp[t] > col) high = t;
83416022c9SBarry Smith         else             low  = t;
8417ab2063SBarry Smith       }
85416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
8617ab2063SBarry Smith         if (rp[i] > col) break;
8717ab2063SBarry Smith         if (rp[i] == col) {
88416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
8917ab2063SBarry Smith           else                  ap[i] = value;
9017ab2063SBarry Smith           goto noinsert;
9117ab2063SBarry Smith         }
9217ab2063SBarry Smith       }
9317ab2063SBarry Smith       if (nonew) goto noinsert;
9417ab2063SBarry Smith       if (nrow >= rmax) {
9517ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
96416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
9717ab2063SBarry Smith         Scalar *new_a;
9817ab2063SBarry Smith 
9917ab2063SBarry Smith         /* malloc new storage space */
100416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1010452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
10217ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
10317ab2063SBarry Smith         new_i   = new_j + new_nz;
10417ab2063SBarry Smith 
10517ab2063SBarry Smith         /* copy over old data into new slots */
10617ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
107416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
108416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
109416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
110416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
11117ab2063SBarry Smith                                                            len*sizeof(int));
112416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
113416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
11417ab2063SBarry Smith                                                            len*sizeof(Scalar));
11517ab2063SBarry Smith         /* free up old matrix storage */
1160452661fSBarry Smith         PetscFree(a->a);
1170452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
118416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
119416022c9SBarry Smith         a->singlemalloc = 1;
12017ab2063SBarry Smith 
12117ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
122416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
123416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
124416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
12517ab2063SBarry Smith       }
126416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
127416022c9SBarry Smith       /* shift up all the later entries in this row */
128416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
12917ab2063SBarry Smith         rp[ii+1] = rp[ii];
13017ab2063SBarry Smith         ap[ii+1] = ap[ii];
13117ab2063SBarry Smith       }
13217ab2063SBarry Smith       rp[i] = col;
13317ab2063SBarry Smith       ap[i] = value;
13417ab2063SBarry Smith       noinsert:;
135416022c9SBarry Smith       low = i + 1;
13617ab2063SBarry Smith     }
13717ab2063SBarry Smith     ailen[row] = nrow;
13817ab2063SBarry Smith   }
13917ab2063SBarry Smith   return 0;
14017ab2063SBarry Smith }
14117ab2063SBarry Smith 
1427eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1437eb43aa7SLois Curfman McInnes {
1447eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
145*b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1467eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1477eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
1487eb43aa7SLois Curfman McInnes 
1497eb43aa7SLois Curfman McInnes   if (!a->assembled) SETERRQ(1,"MatGetValues_SeqAIJ:Not for unassembled matrix");
1507eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
1517eb43aa7SLois Curfman McInnes     row  = im[k];
1527eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
1537eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
1547eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
1557eb43aa7SLois Curfman McInnes     nrow = ailen[row];
1567eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
1577eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
1587eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
1597eb43aa7SLois Curfman McInnes       col = in[l] - shift;
1607eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
1617eb43aa7SLois Curfman McInnes       while (high-low > 5) {
1627eb43aa7SLois Curfman McInnes         t = (low+high)/2;
1637eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
1647eb43aa7SLois Curfman McInnes         else             low  = t;
1657eb43aa7SLois Curfman McInnes       }
1667eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
1677eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
1687eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
169*b49de8d1SLois Curfman McInnes           *v++ = ap[i];
1707eb43aa7SLois Curfman McInnes           goto finished;
1717eb43aa7SLois Curfman McInnes         }
1727eb43aa7SLois Curfman McInnes       }
173*b49de8d1SLois Curfman McInnes       *v++ = zero;
1747eb43aa7SLois Curfman McInnes       finished:;
1757eb43aa7SLois Curfman McInnes     }
1767eb43aa7SLois Curfman McInnes   }
1777eb43aa7SLois Curfman McInnes   return 0;
1787eb43aa7SLois Curfman McInnes }
1797eb43aa7SLois Curfman McInnes 
18017ab2063SBarry Smith #include "draw.h"
18117ab2063SBarry Smith #include "pinclude/pviewer.h"
182416022c9SBarry Smith #include "sysio.h"
18317ab2063SBarry Smith 
184416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
18517ab2063SBarry Smith {
186416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
187416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
18817ab2063SBarry Smith 
189416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
1900452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
191416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
192416022c9SBarry Smith   col_lens[1] = a->m;
193416022c9SBarry Smith   col_lens[2] = a->n;
194416022c9SBarry Smith   col_lens[3] = a->nz;
195416022c9SBarry Smith 
196416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
197416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
198416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
19917ab2063SBarry Smith   }
200416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
2010452661fSBarry Smith   PetscFree(col_lens);
202416022c9SBarry Smith 
203416022c9SBarry Smith   /* store column indices (zero start index) */
204416022c9SBarry Smith   if (a->indexshift) {
205416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
20617ab2063SBarry Smith   }
207416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
208416022c9SBarry Smith   if (a->indexshift) {
209416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
21017ab2063SBarry Smith   }
211416022c9SBarry Smith 
212416022c9SBarry Smith   /* store nonzero values */
213416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
21417ab2063SBarry Smith   return 0;
21517ab2063SBarry Smith }
216416022c9SBarry Smith 
217416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
218416022c9SBarry Smith {
219416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
220416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
22117ab2063SBarry Smith   FILE        *fd;
22217ab2063SBarry Smith   char        *outputname;
22317ab2063SBarry Smith 
224416022c9SBarry Smith   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
225416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
226416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
22717ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
22808480c60SBarry Smith     /* no need to print additional information */ ;
22917ab2063SBarry Smith   }
23017ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
23117ab2063SBarry Smith     int nz, nzalloc, mem;
232416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
233416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
23417ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
23517ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
23617ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
23717ab2063SBarry Smith 
23817ab2063SBarry Smith     for (i=0; i<m; i++) {
239416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
24017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
241416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
242416022c9SBarry Smith                    imag(a->a[j]));
24317ab2063SBarry Smith #else
244416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
24517ab2063SBarry Smith #endif
24617ab2063SBarry Smith       }
24717ab2063SBarry Smith     }
24817ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
24917ab2063SBarry Smith   }
25017ab2063SBarry Smith   else {
25117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
25217ab2063SBarry Smith       fprintf(fd,"row %d:",i);
253416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
25417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
255416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
256416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
25717ab2063SBarry Smith         }
25817ab2063SBarry Smith         else {
259416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
26017ab2063SBarry Smith         }
26117ab2063SBarry Smith #else
262416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
26317ab2063SBarry Smith #endif
26417ab2063SBarry Smith       }
26517ab2063SBarry Smith       fprintf(fd,"\n");
26617ab2063SBarry Smith     }
26717ab2063SBarry Smith   }
26817ab2063SBarry Smith   fflush(fd);
269416022c9SBarry Smith   return 0;
270416022c9SBarry Smith }
271416022c9SBarry Smith 
272416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
273416022c9SBarry Smith {
274416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
275cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
276cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
277d7e8b826SBarry Smith   Draw     draw = (Draw) viewer;
278cddf8d76SBarry Smith   DrawButton  button;
279cddf8d76SBarry Smith 
280416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
281416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
282416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
283416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
284cddf8d76SBarry Smith   color = DRAW_BLUE;
285416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
286cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
287416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
288cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
289cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
290cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
291cddf8d76SBarry Smith #else
292cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
293cddf8d76SBarry Smith #endif
294cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
295cddf8d76SBarry Smith     }
296cddf8d76SBarry Smith   }
297cddf8d76SBarry Smith   color = DRAW_CYAN;
298cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
299cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
300cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
301cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
302cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
303cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
304cddf8d76SBarry Smith     }
305cddf8d76SBarry Smith   }
306cddf8d76SBarry Smith   color = DRAW_RED;
307cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
308cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
309cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
310cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
311cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
312cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
313cddf8d76SBarry Smith #else
314cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
315cddf8d76SBarry Smith #endif
316cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
317416022c9SBarry Smith     }
318416022c9SBarry Smith   }
319416022c9SBarry Smith   DrawFlush(draw);
320cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
321cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
322cddf8d76SBarry Smith 
323cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
324cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
325cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
326cddf8d76SBarry Smith     DrawClear(draw);
327cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
328cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
329cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
330cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
331cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
332cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
333cddf8d76SBarry Smith     w *= scale; h *= scale;
334cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
335cddf8d76SBarry Smith     color = DRAW_BLUE;
336cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
337cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
338cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
339cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
340cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
341cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
342cddf8d76SBarry Smith #else
343cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
344cddf8d76SBarry Smith #endif
345cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
346cddf8d76SBarry Smith       }
347cddf8d76SBarry Smith     }
348cddf8d76SBarry Smith     color = DRAW_CYAN;
349cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
350cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
351cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
352cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
353cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
354cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
355cddf8d76SBarry Smith       }
356cddf8d76SBarry Smith     }
357cddf8d76SBarry Smith     color = DRAW_RED;
358cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
359cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
360cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
361cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
362cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
363cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
364cddf8d76SBarry Smith #else
365cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
366cddf8d76SBarry Smith #endif
367cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
368cddf8d76SBarry Smith       }
369cddf8d76SBarry Smith     }
370cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
371cddf8d76SBarry Smith   }
372416022c9SBarry Smith   return 0;
373416022c9SBarry Smith }
374416022c9SBarry Smith 
375416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
376416022c9SBarry Smith {
377416022c9SBarry Smith   Mat         A = (Mat) obj;
378416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
379416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
380416022c9SBarry Smith 
381416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
382416022c9SBarry Smith   if (!viewer) {
383416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
384416022c9SBarry Smith   }
385416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
386416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
387416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
388416022c9SBarry Smith     }
389416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
390416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
391416022c9SBarry Smith     }
392416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
393416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
394416022c9SBarry Smith     }
395416022c9SBarry Smith   }
396416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
397416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
398416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
39917ab2063SBarry Smith   }
40017ab2063SBarry Smith   return 0;
40117ab2063SBarry Smith }
40241c01911SSatish Balay int Mat_AIJ_CheckInode(Mat);
403416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
40417ab2063SBarry Smith {
405416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
40641c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
407416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
408416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
40917ab2063SBarry Smith 
41017ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
41117ab2063SBarry Smith 
41217ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
413416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
41417ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
41517ab2063SBarry Smith     if (fshift) {
416416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
41717ab2063SBarry Smith       N = ailen[i];
41817ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
41917ab2063SBarry Smith         ip[j-fshift] = ip[j];
42017ab2063SBarry Smith         ap[j-fshift] = ap[j];
42117ab2063SBarry Smith       }
42217ab2063SBarry Smith     }
42317ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
42417ab2063SBarry Smith   }
42517ab2063SBarry Smith   if (m) {
42617ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
42717ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
42817ab2063SBarry Smith   }
42917ab2063SBarry Smith   /* reset ilen and imax for each row */
43017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
43117ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
43217ab2063SBarry Smith   }
433416022c9SBarry Smith   a->nz = ai[m] + shift;
43417ab2063SBarry Smith 
43517ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
436416022c9SBarry Smith   if (fshift && a->diag) {
4370452661fSBarry Smith     PetscFree(a->diag);
438416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
439416022c9SBarry Smith     a->diag = 0;
44017ab2063SBarry Smith   }
44176dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
44241c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
443416022c9SBarry Smith   a->assembled = 1;
44417ab2063SBarry Smith   return 0;
44517ab2063SBarry Smith }
44617ab2063SBarry Smith 
447416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
44817ab2063SBarry Smith {
449416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
450cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
45117ab2063SBarry Smith   return 0;
45217ab2063SBarry Smith }
453416022c9SBarry Smith 
45417ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
45517ab2063SBarry Smith {
456416022c9SBarry Smith   Mat        A  = (Mat) obj;
457416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
458d5d45c9bSBarry Smith 
45917ab2063SBarry Smith #if defined(PETSC_LOG)
460416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
46117ab2063SBarry Smith #endif
4620452661fSBarry Smith   PetscFree(a->a);
4630452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4640452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
4650452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4660452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
4670452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
46876dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
4690452661fSBarry Smith   PetscFree(a);
470416022c9SBarry Smith   PLogObjectDestroy(A);
4710452661fSBarry Smith   PetscHeaderDestroy(A);
47217ab2063SBarry Smith   return 0;
47317ab2063SBarry Smith }
47417ab2063SBarry Smith 
475416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
47617ab2063SBarry Smith {
47717ab2063SBarry Smith   return 0;
47817ab2063SBarry Smith }
47917ab2063SBarry Smith 
480416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
48117ab2063SBarry Smith {
482416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
483416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
484416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
485416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
486416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
487e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
488e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
489e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
490e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
491df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
492df719601SLois Curfman McInnes   else if (op == COLUMN_ORIENTED)
4934d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");}
494df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
4954d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
496e2f28af5SBarry Smith   else
4974d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
49817ab2063SBarry Smith   return 0;
49917ab2063SBarry Smith }
50017ab2063SBarry Smith 
501416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
50217ab2063SBarry Smith {
503416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
504416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
50517ab2063SBarry Smith   Scalar     *x, zero = 0.0;
50617ab2063SBarry Smith 
507416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
50817ab2063SBarry Smith   VecSet(&zero,v);
50917ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
510416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
511416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
512416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
513416022c9SBarry Smith       if (a->j[j]+shift == i) {
514416022c9SBarry Smith         x[i] = a->a[j];
51517ab2063SBarry Smith         break;
51617ab2063SBarry Smith       }
51717ab2063SBarry Smith     }
51817ab2063SBarry Smith   }
51917ab2063SBarry Smith   return 0;
52017ab2063SBarry Smith }
52117ab2063SBarry Smith 
52217ab2063SBarry Smith /* -------------------------------------------------------*/
52317ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
52417ab2063SBarry Smith /* -------------------------------------------------------*/
525416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
52617ab2063SBarry Smith {
527416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
52817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
529416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
53017ab2063SBarry Smith 
531416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
53217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
533cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
53417ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
53517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
536416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
537416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
538416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
53917ab2063SBarry Smith     alpha = x[i];
54017ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
54117ab2063SBarry Smith   }
542416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
54317ab2063SBarry Smith   return 0;
54417ab2063SBarry Smith }
545d5d45c9bSBarry Smith 
546416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
54717ab2063SBarry Smith {
548416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
54917ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
550416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
55117ab2063SBarry Smith 
552416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
55317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
55417ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
55517ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
55617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
557416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
558416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
559416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
56017ab2063SBarry Smith     alpha = x[i];
56117ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
56217ab2063SBarry Smith   }
56317ab2063SBarry Smith   return 0;
56417ab2063SBarry Smith }
56517ab2063SBarry Smith 
566416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
56717ab2063SBarry Smith {
568416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
56917ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
570416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
57117ab2063SBarry Smith 
572416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
57317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
57417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
575416022c9SBarry Smith   idx  = a->j;
576416022c9SBarry Smith   v    = a->a;
577416022c9SBarry Smith   ii   = a->i;
57817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
579416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
58017ab2063SBarry Smith     sum  = 0.0;
58117ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
58217ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
583416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
58417ab2063SBarry Smith     y[i] = sum;
58517ab2063SBarry Smith   }
586416022c9SBarry Smith   PLogFlops(2*a->nz - m);
58717ab2063SBarry Smith   return 0;
58817ab2063SBarry Smith }
58917ab2063SBarry Smith 
590416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
59117ab2063SBarry Smith {
592416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
59317ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
594cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
59517ab2063SBarry Smith 
59648d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
59717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
59817ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
599cddf8d76SBarry Smith   idx  = a->j;
600cddf8d76SBarry Smith   v    = a->a;
601cddf8d76SBarry Smith   ii   = a->i;
60217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
603cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
60417ab2063SBarry Smith     sum  = y[i];
605cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
606cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
60717ab2063SBarry Smith     z[i] = sum;
60817ab2063SBarry Smith   }
609416022c9SBarry Smith   PLogFlops(2*a->nz);
61017ab2063SBarry Smith   return 0;
61117ab2063SBarry Smith }
61217ab2063SBarry Smith 
61317ab2063SBarry Smith /*
61417ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
61517ab2063SBarry Smith */
61617ab2063SBarry Smith 
617416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
61817ab2063SBarry Smith {
619416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
620416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
62117ab2063SBarry Smith 
622416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
6230452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
624416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
625416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
626416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
627416022c9SBarry Smith       if (a->j[j]+shift == i) {
62817ab2063SBarry Smith         diag[i] = j - shift;
62917ab2063SBarry Smith         break;
63017ab2063SBarry Smith       }
63117ab2063SBarry Smith     }
63217ab2063SBarry Smith   }
633416022c9SBarry Smith   a->diag = diag;
63417ab2063SBarry Smith   return 0;
63517ab2063SBarry Smith }
63617ab2063SBarry Smith 
637416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
63817ab2063SBarry Smith                            double fshift,int its,Vec xx)
63917ab2063SBarry Smith {
640416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
641416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
642d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
64317ab2063SBarry Smith 
64417ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
645416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
646416022c9SBarry Smith   diag = a->diag;
647416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
64817ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
64917ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
65017ab2063SBarry Smith     bs = b + shift;
65117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
652416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
653416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
654416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
655416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
65617ab2063SBarry Smith         sum  = b[i]*d/omega;
65717ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
65817ab2063SBarry Smith         x[i] = sum;
65917ab2063SBarry Smith     }
66017ab2063SBarry Smith     return 0;
66117ab2063SBarry Smith   }
66217ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
66317ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
66417ab2063SBarry Smith   }
665416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
66617ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
66717ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
66817ab2063SBarry Smith 
66917ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
67017ab2063SBarry Smith 
67117ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
67217ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
67317ab2063SBarry Smith     is the relaxation factor.
67417ab2063SBarry Smith     */
6750452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
67617ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
67717ab2063SBarry Smith 
67817ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
67917ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
680416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
681416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
682416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
683416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
68417ab2063SBarry Smith       sum  = b[i];
68517ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
68617ab2063SBarry Smith       x[i] = omega*(sum/d);
68717ab2063SBarry Smith     }
68817ab2063SBarry Smith 
68917ab2063SBarry Smith     /*  t = b - (2*E - D)x */
690416022c9SBarry Smith     v = a->a;
69117ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
69217ab2063SBarry Smith 
69317ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
694416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
695416022c9SBarry Smith     diag = a->diag;
69617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
697416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
698416022c9SBarry Smith       n    = diag[i] - a->i[i];
699416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
700416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
70117ab2063SBarry Smith       sum  = t[i];
70217ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
70317ab2063SBarry Smith       t[i] = omega*(sum/d);
70417ab2063SBarry Smith     }
70517ab2063SBarry Smith 
70617ab2063SBarry Smith     /*  x = x + t */
70717ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7080452661fSBarry Smith     PetscFree(t);
70917ab2063SBarry Smith     return 0;
71017ab2063SBarry Smith   }
71117ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
71217ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
71317ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
714416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
715416022c9SBarry Smith         n    = diag[i] - a->i[i];
716416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
717416022c9SBarry Smith         v    = a->a + a->i[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       xb = x;
72317ab2063SBarry Smith     }
72417ab2063SBarry Smith     else xb = b;
72517ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
72617ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
72717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
728416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
72917ab2063SBarry Smith       }
73017ab2063SBarry Smith     }
73117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
73217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
733416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
734416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
735416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
736416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
73717ab2063SBarry Smith         sum  = xb[i];
73817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
73917ab2063SBarry Smith         x[i] = omega*(sum/d);
74017ab2063SBarry Smith       }
74117ab2063SBarry Smith     }
74217ab2063SBarry Smith     its--;
74317ab2063SBarry Smith   }
74417ab2063SBarry Smith   while (its--) {
74517ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
74617ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
747416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
748416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
749416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
750416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
75117ab2063SBarry Smith         sum  = b[i];
75217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
75317ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
75417ab2063SBarry Smith       }
75517ab2063SBarry Smith     }
75617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
75717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
758416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
759416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
760416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
761416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
76217ab2063SBarry Smith         sum  = b[i];
76317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
76417ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
76517ab2063SBarry Smith       }
76617ab2063SBarry Smith     }
76717ab2063SBarry Smith   }
76817ab2063SBarry Smith   return 0;
76917ab2063SBarry Smith }
77017ab2063SBarry Smith 
771d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
77217ab2063SBarry Smith {
773416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
774416022c9SBarry Smith   *nz      = a->nz;
775416022c9SBarry Smith   *nzalloc = a->maxnz;
776416022c9SBarry Smith   *mem     = (int)A->mem;
77717ab2063SBarry Smith   return 0;
77817ab2063SBarry Smith }
77917ab2063SBarry Smith 
78017ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
78117ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
78217ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
78317ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
78417ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
78517ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
78617ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
78717ab2063SBarry Smith 
78817ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
78917ab2063SBarry Smith {
790416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
791416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
79217ab2063SBarry Smith 
79317ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
79417ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
79517ab2063SBarry Smith   if (diag) {
79617ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
797416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
798416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
799416022c9SBarry Smith         a->ilen[rows[i]] = 1;
800416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
801416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
80217ab2063SBarry Smith       }
80317ab2063SBarry Smith       else {
80417ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
80517ab2063SBarry Smith         CHKERRQ(ierr);
80617ab2063SBarry Smith       }
80717ab2063SBarry Smith     }
80817ab2063SBarry Smith   }
80917ab2063SBarry Smith   else {
81017ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
811416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
812416022c9SBarry Smith       a->ilen[rows[i]] = 0;
81317ab2063SBarry Smith     }
81417ab2063SBarry Smith   }
81517ab2063SBarry Smith   ISRestoreIndices(is,&rows);
81617ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
81717ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
81817ab2063SBarry Smith   return 0;
81917ab2063SBarry Smith }
82017ab2063SBarry Smith 
821416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
82217ab2063SBarry Smith {
823416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
824416022c9SBarry Smith   *m = a->m; *n = a->n;
82517ab2063SBarry Smith   return 0;
82617ab2063SBarry Smith }
82717ab2063SBarry Smith 
828416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
82917ab2063SBarry Smith {
830416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
831416022c9SBarry Smith   *m = 0; *n = a->m;
83217ab2063SBarry Smith   return 0;
83317ab2063SBarry Smith }
834416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
83517ab2063SBarry Smith {
836416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
837416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
83817ab2063SBarry Smith 
839416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
84017ab2063SBarry Smith 
841416022c9SBarry Smith   if (!a->assembled) {
842416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
843416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
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 
871416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
87217ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
873416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
87417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
87517ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
87617ab2063SBarry Smith #else
87717ab2063SBarry Smith       sum += (*v)*(*v); v++;
87817ab2063SBarry Smith #endif
87917ab2063SBarry Smith     }
88017ab2063SBarry Smith     *norm = sqrt(sum);
88117ab2063SBarry Smith   }
88217ab2063SBarry Smith   else if (type == NORM_1) {
88317ab2063SBarry Smith     double *tmp;
884416022c9SBarry Smith     int    *jj = a->j;
8850452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
886cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
88717ab2063SBarry Smith     *norm = 0.0;
888416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
889a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
89017ab2063SBarry Smith     }
891416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
89217ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
89317ab2063SBarry Smith     }
8940452661fSBarry Smith     PetscFree(tmp);
89517ab2063SBarry Smith   }
89617ab2063SBarry Smith   else if (type == NORM_INFINITY) {
89717ab2063SBarry Smith     *norm = 0.0;
898416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
899416022c9SBarry Smith       v = a->a + a->i[j] + shift;
90017ab2063SBarry Smith       sum = 0.0;
901416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
902cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
90317ab2063SBarry Smith       }
90417ab2063SBarry Smith       if (sum > *norm) *norm = sum;
90517ab2063SBarry Smith     }
90617ab2063SBarry Smith   }
90717ab2063SBarry Smith   else {
90848d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
90917ab2063SBarry Smith   }
91017ab2063SBarry Smith   return 0;
91117ab2063SBarry Smith }
91217ab2063SBarry Smith 
913416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
91417ab2063SBarry Smith {
915416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
916416022c9SBarry Smith   Mat        C;
917416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
918416022c9SBarry Smith   int        shift = a->indexshift;
919d5d45c9bSBarry Smith   Scalar     *array = a->a;
92017ab2063SBarry Smith 
921416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat 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 
942416022c9SBarry Smith   if (B) {
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 
959416022c9SBarry Smith static int MatScale_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 
96548d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
96617ab2063SBarry Smith   if (ll) {
96717ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
968416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
969416022c9SBarry Smith     v = a->a;
97017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
97117ab2063SBarry Smith       x = l[i];
972416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
97317ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
97417ab2063SBarry Smith     }
97517ab2063SBarry Smith   }
97617ab2063SBarry Smith   if (rr) {
97717ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
978416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
979416022c9SBarry Smith     v = a->a; jj = a->j;
98017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
98117ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
98217ab2063SBarry Smith     }
98317ab2063SBarry Smith   }
98417ab2063SBarry Smith   return 0;
98517ab2063SBarry Smith }
98617ab2063SBarry Smith 
987cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
98817ab2063SBarry Smith {
989db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
99002834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
991a2744918SBarry Smith   register int sum,lensi;
99202834360SBarry Smith   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
993a2744918SBarry Smith   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
994db02288aSLois Curfman McInnes   Scalar       *vwork,*a_new;
995416022c9SBarry Smith   Mat          C;
99617ab2063SBarry Smith 
997416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
99817ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
99917ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
100017ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
100117ab2063SBarry Smith 
100202834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
100302834360SBarry Smith     /* special case of contiguous rows */
10040452661fSBarry Smith     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
100502834360SBarry Smith     starts = lens + ncols;
100602834360SBarry Smith     /* loop over new rows determining lens and starting points */
100702834360SBarry Smith     for (i=0; i<nrows; i++) {
1008a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1009a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
101002834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1011a2744918SBarry Smith         if (aj[k] >= first) {
101202834360SBarry Smith           starts[i] = k;
101302834360SBarry Smith           break;
101402834360SBarry Smith 	}
101502834360SBarry Smith       }
1016a2744918SBarry Smith       sum = 0;
101702834360SBarry Smith       while (k < kend) {
1018a2744918SBarry Smith         if (aj[k++] >= first+ncols) break;
1019a2744918SBarry Smith         sum++;
102002834360SBarry Smith       }
1021a2744918SBarry Smith       lens[i] = sum;
102202834360SBarry Smith     }
102302834360SBarry Smith     /* create submatrix */
1024cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
102508480c60SBarry Smith       int n_cols,n_rows;
102608480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
102708480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
102808480c60SBarry Smith       MatZeroEntries(*B);
102908480c60SBarry Smith       C = *B;
103008480c60SBarry Smith     }
103108480c60SBarry Smith     else {
103202834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
103308480c60SBarry Smith     }
1034db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1035db02288aSLois Curfman McInnes 
103602834360SBarry Smith     /* loop over rows inserting into submatrix */
1037db02288aSLois Curfman McInnes     a_new    = c->a;
1038db02288aSLois Curfman McInnes     j_new    = c->j;
1039db02288aSLois Curfman McInnes     i_new    = c->i;
1040db02288aSLois Curfman McInnes     i_new[0] = -shift;
104102834360SBarry Smith     for (i=0; i<nrows; i++) {
1042a2744918SBarry Smith       ii    = starts[i];
1043a2744918SBarry Smith       lensi = lens[i];
1044a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1045a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
104602834360SBarry Smith       }
1047a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1048a2744918SBarry Smith       a_new      += lensi;
1049a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1050a2744918SBarry Smith       c->ilen[i]  = lensi;
105102834360SBarry Smith     }
10520452661fSBarry Smith     PetscFree(lens);
105302834360SBarry Smith   }
105402834360SBarry Smith   else {
105502834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
10560452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
105702834360SBarry Smith     ssmap = smap + shift;
10587eb43aa7SLois Curfman McInnes     cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork);
105902834360SBarry Smith     lens  = cwork + ncols;
10600452661fSBarry Smith     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1061cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
106217ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
106302834360SBarry Smith     /* determine lens of each row */
106402834360SBarry Smith     for (i=0; i<nrows; i++) {
106502834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
106602834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
106702834360SBarry Smith       lens[i] = 0;
106802834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
106902834360SBarry Smith         if (ssmap[a->j[k]]) {
107002834360SBarry Smith           lens[i]++;
107102834360SBarry Smith         }
107202834360SBarry Smith       }
107302834360SBarry Smith     }
107417ab2063SBarry Smith     /* Create and fill new matrix */
1075a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
107608480c60SBarry Smith       int n_cols,n_rows;
107708480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
107808480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
107908480c60SBarry Smith       MatZeroEntries(*B);
108008480c60SBarry Smith       C = *B;
108108480c60SBarry Smith     }
108208480c60SBarry Smith     else {
108302834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
108408480c60SBarry Smith     }
108517ab2063SBarry Smith     for (i=0; i<nrows; i++) {
108617ab2063SBarry Smith       nznew  = 0;
1087416022c9SBarry Smith       kstart = a->i[irow[i]]+shift;
1088416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
108917ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
109002834360SBarry Smith         if (ssmap[a->j[k]]) {
109102834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
1092416022c9SBarry Smith           vwork[nznew++] = a->a[k];
109317ab2063SBarry Smith         }
109417ab2063SBarry Smith       }
1095416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
109617ab2063SBarry Smith     }
109702834360SBarry Smith     /* Free work space */
109802834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
10990452661fSBarry Smith     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
110002834360SBarry Smith   }
1101416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1102416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
110317ab2063SBarry Smith 
110417ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1105416022c9SBarry Smith   *B = C;
110617ab2063SBarry Smith   return 0;
110717ab2063SBarry Smith }
110817ab2063SBarry Smith 
1109a871dcd8SBarry Smith /*
111063b91edcSBarry Smith      note: This can only work for identity for row and col. It would
111163b91edcSBarry Smith    be good to check this and otherwise generate an error.
1112a871dcd8SBarry Smith */
111363b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1114a871dcd8SBarry Smith {
111563b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
111608480c60SBarry Smith   int        ierr;
111763b91edcSBarry Smith   Mat        outA;
111863b91edcSBarry Smith 
1119a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1120a871dcd8SBarry Smith 
112163b91edcSBarry Smith   outA          = inA;
112263b91edcSBarry Smith   inA->factor   = FACTOR_LU;
112363b91edcSBarry Smith   a->row        = row;
112463b91edcSBarry Smith   a->col        = col;
112563b91edcSBarry Smith 
11260452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
112763b91edcSBarry Smith 
112808480c60SBarry Smith   if (!a->diag) {
112908480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
113063b91edcSBarry Smith   }
113163b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1132a871dcd8SBarry Smith   return 0;
1133a871dcd8SBarry Smith }
1134a871dcd8SBarry Smith 
1135cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1136cddf8d76SBarry Smith                                     Mat **B)
1137cddf8d76SBarry Smith {
1138cddf8d76SBarry Smith   int ierr,i;
1139cddf8d76SBarry Smith 
1140cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
11410452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1142cddf8d76SBarry Smith   }
1143cddf8d76SBarry Smith 
1144cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1145cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1146cddf8d76SBarry Smith   }
1147cddf8d76SBarry Smith   return 0;
1148cddf8d76SBarry Smith }
1149cddf8d76SBarry Smith 
11504dcbc457SBarry Smith static int MatIncreaseOverlap_SeqAIJ(Mat A, int n, IS *is, int ov)
11514dcbc457SBarry Smith {
11524dcbc457SBarry Smith   int i,m,*idx,ierr;
11534dcbc457SBarry Smith 
11544dcbc457SBarry Smith   for ( i=0; i<n; i++ ) {
11554dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr);
11564dcbc457SBarry Smith     ISGetLocalSize(is[i],&m);
11574dcbc457SBarry Smith   }
11584dcbc457SBarry Smith   SETERRQ(1,"MatIncreaseOverlap_SeqAIJ:Not implemented");
11594dcbc457SBarry Smith }
116017ab2063SBarry Smith /* -------------------------------------------------------------------*/
116117ab2063SBarry Smith 
116217ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
116317ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1164416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1165416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
116617ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
116717ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
116817ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
116917ab2063SBarry Smith        MatRelax_SeqAIJ,
117017ab2063SBarry Smith        MatTranspose_SeqAIJ,
117117ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
117217ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
117317ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
117417ab2063SBarry Smith        MatCompress_SeqAIJ,
117517ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
117617ab2063SBarry Smith        MatGetReordering_SeqAIJ,
117717ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
117817ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
117917ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
118017ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1181416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
1182a871dcd8SBarry Smith        MatCopyPrivate_SeqAIJ,0,0,
1183cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
11847eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
11857eb43aa7SLois Curfman McInnes        MatGetValues_SeqAIJ};
118617ab2063SBarry Smith 
118717ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
118817ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
118917ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
119017ab2063SBarry Smith 
119117ab2063SBarry Smith /*@C
119217ab2063SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
119317ab2063SBarry Smith    (the default uniprocessor PETSc format).
119417ab2063SBarry Smith 
119517ab2063SBarry Smith    Input Parameters:
119617ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
119717ab2063SBarry Smith .  m - number of rows
119817ab2063SBarry Smith .  n - number of columns
119917ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
120017ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
120117ab2063SBarry Smith 
120217ab2063SBarry Smith    Output Parameter:
1203416022c9SBarry Smith .  A - the matrix
120417ab2063SBarry Smith 
120517ab2063SBarry Smith    Notes:
120617ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
120717ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
12080002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
12090002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
121017ab2063SBarry Smith 
121117ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
121217ab2063SBarry Smith    Set both nz and nnz to zero for PETSc to control dynamic memory
121317ab2063SBarry Smith    allocation.
121417ab2063SBarry Smith 
121517ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse
121617ab2063SBarry Smith 
121717ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
121817ab2063SBarry Smith @*/
1219416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
122017ab2063SBarry Smith {
1221416022c9SBarry Smith   Mat        B;
1222416022c9SBarry Smith   Mat_SeqAIJ *b;
122317ab2063SBarry Smith   int        i,len,ierr;
1224d5d45c9bSBarry Smith 
1225416022c9SBarry Smith   *A      = 0;
12260452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1227416022c9SBarry Smith   PLogObjectCreate(B);
12280452661fSBarry Smith   B->data               = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1229416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1230416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1231416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1232416022c9SBarry Smith   B->factor           = 0;
1233416022c9SBarry Smith   B->lupivotthreshold = 1.0;
1234416022c9SBarry Smith   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1235416022c9SBarry Smith   b->row              = 0;
1236416022c9SBarry Smith   b->col              = 0;
1237416022c9SBarry Smith   b->indexshift       = 0;
1238416022c9SBarry Smith   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
123917ab2063SBarry Smith 
1240416022c9SBarry Smith   b->m       = m;
1241416022c9SBarry Smith   b->n       = n;
12420452661fSBarry Smith   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
124317ab2063SBarry Smith   if (!nnz) {
124417ab2063SBarry Smith     if (nz <= 0) nz = 1;
1245416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
124617ab2063SBarry Smith     nz = nz*m;
124717ab2063SBarry Smith   }
124817ab2063SBarry Smith   else {
124917ab2063SBarry Smith     nz = 0;
1250416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
125117ab2063SBarry Smith   }
125217ab2063SBarry Smith 
125317ab2063SBarry Smith   /* allocate the matrix space */
1254416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
12550452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1256416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1257cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1258416022c9SBarry Smith   b->i  = b->j + nz;
1259416022c9SBarry Smith   b->singlemalloc = 1;
126017ab2063SBarry Smith 
1261416022c9SBarry Smith   b->i[0] = -b->indexshift;
126217ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1263416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
126417ab2063SBarry Smith   }
126517ab2063SBarry Smith 
1266416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
12670452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1268416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1269416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
127017ab2063SBarry Smith 
1271416022c9SBarry Smith   b->nz               = 0;
1272416022c9SBarry Smith   b->maxnz            = nz;
1273416022c9SBarry Smith   b->sorted           = 0;
1274416022c9SBarry Smith   b->roworiented      = 1;
1275416022c9SBarry Smith   b->nonew            = 0;
1276416022c9SBarry Smith   b->diag             = 0;
1277416022c9SBarry Smith   b->assembled        = 0;
1278416022c9SBarry Smith   b->solve_work       = 0;
127971bd300dSLois Curfman McInnes   b->spptr            = 0;
1280754ec7b1SSatish Balay   b->inode.node_count = 0;
1281754ec7b1SSatish Balay   b->inode.size       = 0;
128217ab2063SBarry Smith 
1283416022c9SBarry Smith   *A = B;
128417ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_superlu")) {
1285416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
128617ab2063SBarry Smith   }
128717ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_essl")) {
1288416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
128917ab2063SBarry Smith   }
129017ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_dxml")) {
1291416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1292416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
129317ab2063SBarry Smith   }
129417ab2063SBarry Smith 
129517ab2063SBarry Smith   return 0;
129617ab2063SBarry Smith }
129717ab2063SBarry Smith 
129808480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues)
129917ab2063SBarry Smith {
1300416022c9SBarry Smith   Mat        C;
1301416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
130208480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
130317ab2063SBarry Smith 
13044043dd9cSLois Curfman McInnes   *B = 0;
1305416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
13060452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1307416022c9SBarry Smith   PLogObjectCreate(C);
13080452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
130941c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1310416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1311416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1312416022c9SBarry Smith   C->factor     = A->factor;
1313416022c9SBarry Smith   c->row        = 0;
1314416022c9SBarry Smith   c->col        = 0;
1315416022c9SBarry Smith   c->indexshift = shift;
131617ab2063SBarry Smith 
1317416022c9SBarry Smith   c->m          = a->m;
1318416022c9SBarry Smith   c->n          = a->n;
131917ab2063SBarry Smith 
13200452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
13210452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
132217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1323416022c9SBarry Smith     c->imax[i] = a->imax[i];
1324416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
132517ab2063SBarry Smith   }
132617ab2063SBarry Smith 
132717ab2063SBarry Smith   /* allocate the matrix space */
1328416022c9SBarry Smith   c->singlemalloc = 1;
1329416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
13300452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1331416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1332416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1333416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
133417ab2063SBarry Smith   if (m > 0) {
1335416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
133608480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1337416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
133817ab2063SBarry Smith     }
133908480c60SBarry Smith   }
134017ab2063SBarry Smith 
1341416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1342416022c9SBarry Smith   c->sorted      = a->sorted;
1343416022c9SBarry Smith   c->roworiented = a->roworiented;
1344416022c9SBarry Smith   c->nonew       = a->nonew;
134517ab2063SBarry Smith 
1346416022c9SBarry Smith   if (a->diag) {
13470452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1348416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
134917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1350416022c9SBarry Smith       c->diag[i] = a->diag[i];
135117ab2063SBarry Smith     }
135217ab2063SBarry Smith   }
1353416022c9SBarry Smith   else c->diag        = 0;
1354754ec7b1SSatish Balay   if( a->inode.size){
1355754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1356754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1357754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1358754ec7b1SSatish Balay   } else {
1359754ec7b1SSatish Balay     c->inode.size       = 0;
1360754ec7b1SSatish Balay     c->inode.node_count = 0;
1361754ec7b1SSatish Balay   }
1362416022c9SBarry Smith   c->assembled        = 1;
1363416022c9SBarry Smith   c->nz               = a->nz;
1364416022c9SBarry Smith   c->maxnz            = a->maxnz;
1365416022c9SBarry Smith   c->solve_work       = 0;
136676dd722bSSatish Balay   c->spptr            = 0;      /* Dangerous -I'm throwing away a->spptr */
1367754ec7b1SSatish Balay 
1368416022c9SBarry Smith   *B = C;
136917ab2063SBarry Smith   return 0;
137017ab2063SBarry Smith }
137117ab2063SBarry Smith 
1372416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
137317ab2063SBarry Smith {
1374416022c9SBarry Smith   Mat_SeqAIJ   *a;
1375416022c9SBarry Smith   Mat          B;
137617699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
137717ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
137817ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
137917ab2063SBarry Smith 
138017699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
138117699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
138217ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1383416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
138448d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
138517ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
138617ab2063SBarry Smith 
138717ab2063SBarry Smith   /* read in row lengths */
13880452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1389416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
139017ab2063SBarry Smith 
139117ab2063SBarry Smith   /* create our matrix */
1392416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1393416022c9SBarry Smith   B = *A;
1394416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1395416022c9SBarry Smith   shift = a->indexshift;
139617ab2063SBarry Smith 
139717ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1398416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
139917ab2063SBarry Smith   if (shift) {
140017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1401416022c9SBarry Smith       a->j[i] += 1;
140217ab2063SBarry Smith     }
140317ab2063SBarry Smith   }
140417ab2063SBarry Smith 
140517ab2063SBarry Smith   /* read in nonzero values */
1406416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
140717ab2063SBarry Smith 
140817ab2063SBarry Smith   /* set matrix "i" values */
1409416022c9SBarry Smith   a->i[0] = -shift;
141017ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1411416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1412416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
141317ab2063SBarry Smith   }
14140452661fSBarry Smith   PetscFree(rowlengths);
141517ab2063SBarry Smith 
1416416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1417416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
141817ab2063SBarry Smith   return 0;
141917ab2063SBarry Smith }
142017ab2063SBarry Smith 
142117ab2063SBarry Smith 
142217ab2063SBarry Smith 
1423