xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 08480c60afa5ef1d2e4e27b9ebdf48b02c6a2186)
117ab2063SBarry Smith #ifndef lint
2*08480c60SBarry Smith static char vcid[] = "$Id: aij.c,v 1.105 1995/10/20 02:00:02 curfman Exp bsmith $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
517ab2063SBarry Smith #include "aij.h"
617ab2063SBarry Smith #include "vec/vecimpl.h"
717ab2063SBarry Smith #include "inline/spops.h"
817ab2063SBarry Smith 
917ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**);
1017ab2063SBarry Smith 
11416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
1217ab2063SBarry Smith {
13416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1417ab2063SBarry Smith   int        ierr, *ia, *ja;
1517ab2063SBarry Smith 
16416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix");
1717ab2063SBarry Smith 
18416022c9SBarry Smith   ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr);
19416022c9SBarry Smith   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
2017ab2063SBarry Smith   PETSCFREE(ia); PETSCFREE(ja);
2117ab2063SBarry Smith   return 0;
2217ab2063SBarry Smith }
2317ab2063SBarry Smith 
24416022c9SBarry Smith #define CHUNKSIZE   10
2517ab2063SBarry Smith 
2617ab2063SBarry Smith /* This version has row oriented v  */
27416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
2817ab2063SBarry Smith {
29416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
30416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
31416022c9SBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen;
32416022c9SBarry Smith   int        *aj = a->j, nonew = a->nonew;
33416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
34416022c9SBarry Smith   int        shift = a->indexshift;
3517ab2063SBarry Smith 
3617ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
37416022c9SBarry Smith     row  = im[k];
3817ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
39416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
4017ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
4117ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
42416022c9SBarry Smith     low = 0;
4317ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
44416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
45416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
46416022c9SBarry Smith       col = in[l] - shift; value = *v++;
47416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
48416022c9SBarry Smith       while (high-low > 5) {
49416022c9SBarry Smith         t = (low+high)/2;
50416022c9SBarry Smith         if (rp[t] > col) high = t;
51416022c9SBarry Smith         else             low  = t;
5217ab2063SBarry Smith       }
53416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
5417ab2063SBarry Smith         if (rp[i] > col) break;
5517ab2063SBarry Smith         if (rp[i] == col) {
56416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
5717ab2063SBarry Smith           else                  ap[i] = value;
5817ab2063SBarry Smith           goto noinsert;
5917ab2063SBarry Smith         }
6017ab2063SBarry Smith       }
6117ab2063SBarry Smith       if (nonew) goto noinsert;
6217ab2063SBarry Smith       if (nrow >= rmax) {
6317ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
64416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
6517ab2063SBarry Smith         Scalar *new_a;
6617ab2063SBarry Smith 
6717ab2063SBarry Smith         /* malloc new storage space */
68416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
6917ab2063SBarry Smith         new_a   = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(new_a);
7017ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
7117ab2063SBarry Smith         new_i   = new_j + new_nz;
7217ab2063SBarry Smith 
7317ab2063SBarry Smith         /* copy over old data into new slots */
7417ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
75416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
76416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
77416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
78416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
7917ab2063SBarry Smith                                                            len*sizeof(int));
80416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
81416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
8217ab2063SBarry Smith                                                            len*sizeof(Scalar));
8317ab2063SBarry Smith         /* free up old matrix storage */
84416022c9SBarry Smith         PETSCFREE(a->a);
85416022c9SBarry Smith         if (!a->singlemalloc) {PETSCFREE(a->i);PETSCFREE(a->j);}
86416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
87416022c9SBarry Smith         a->singlemalloc = 1;
8817ab2063SBarry Smith 
8917ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
90416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
91416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
92416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
9317ab2063SBarry Smith       }
94416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
95416022c9SBarry Smith       /* shift up all the later entries in this row */
96416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
9717ab2063SBarry Smith         rp[ii+1] = rp[ii];
9817ab2063SBarry Smith         ap[ii+1] = ap[ii];
9917ab2063SBarry Smith       }
10017ab2063SBarry Smith       rp[i] = col;
10117ab2063SBarry Smith       ap[i] = value;
10217ab2063SBarry Smith       noinsert:;
103416022c9SBarry Smith       low = i + 1;
10417ab2063SBarry Smith     }
10517ab2063SBarry Smith     ailen[row] = nrow;
10617ab2063SBarry Smith   }
10717ab2063SBarry Smith   return 0;
10817ab2063SBarry Smith }
10917ab2063SBarry Smith 
11017ab2063SBarry Smith #include "draw.h"
11117ab2063SBarry Smith #include "pinclude/pviewer.h"
112416022c9SBarry Smith #include "sysio.h"
11317ab2063SBarry Smith 
114416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
11517ab2063SBarry Smith {
116416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
117416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
11817ab2063SBarry Smith 
119416022c9SBarry Smith   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
120a871dcd8SBarry Smith   col_lens = (int *) PETSCMALLOC( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
121416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
122416022c9SBarry Smith   col_lens[1] = a->m;
123416022c9SBarry Smith   col_lens[2] = a->n;
124416022c9SBarry Smith   col_lens[3] = a->nz;
125416022c9SBarry Smith 
126416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
127416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
128416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
12917ab2063SBarry Smith   }
130416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
131416022c9SBarry Smith   PETSCFREE(col_lens);
132416022c9SBarry Smith 
133416022c9SBarry Smith   /* store column indices (zero start index) */
134416022c9SBarry Smith   if (a->indexshift) {
135416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
13617ab2063SBarry Smith   }
137416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
138416022c9SBarry Smith   if (a->indexshift) {
139416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
14017ab2063SBarry Smith   }
141416022c9SBarry Smith 
142416022c9SBarry Smith   /* store nonzero values */
143416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
14417ab2063SBarry Smith   return 0;
14517ab2063SBarry Smith }
146416022c9SBarry Smith 
147416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
148416022c9SBarry Smith {
149416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
150416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
15117ab2063SBarry Smith   FILE        *fd;
15217ab2063SBarry Smith   char        *outputname;
15317ab2063SBarry Smith 
154416022c9SBarry Smith   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
155416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
156416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
15717ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
158*08480c60SBarry Smith     /* no need to print additional information */ ;
15917ab2063SBarry Smith   }
16017ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
16117ab2063SBarry Smith     int nz, nzalloc, mem;
162416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
163416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
16417ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
16517ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
16617ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
16717ab2063SBarry Smith 
16817ab2063SBarry Smith     for (i=0; i<m; i++) {
169416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
17017ab2063SBarry Smith #if defined(PETSC_COMPLEX)
171416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
172416022c9SBarry Smith                    imag(a->a[j]));
17317ab2063SBarry Smith #else
174416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
17517ab2063SBarry Smith #endif
17617ab2063SBarry Smith       }
17717ab2063SBarry Smith     }
17817ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
17917ab2063SBarry Smith   }
18017ab2063SBarry Smith   else {
18117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
18217ab2063SBarry Smith       fprintf(fd,"row %d:",i);
183416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
18417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
185416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
186416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
18717ab2063SBarry Smith         }
18817ab2063SBarry Smith         else {
189416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
19017ab2063SBarry Smith         }
19117ab2063SBarry Smith #else
192416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
19317ab2063SBarry Smith #endif
19417ab2063SBarry Smith       }
19517ab2063SBarry Smith       fprintf(fd,"\n");
19617ab2063SBarry Smith     }
19717ab2063SBarry Smith   }
19817ab2063SBarry Smith   fflush(fd);
199416022c9SBarry Smith   return 0;
200416022c9SBarry Smith }
201416022c9SBarry Smith 
202416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
203416022c9SBarry Smith {
204416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
205416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift;
206416022c9SBarry Smith   double      xl,yl,xr,yr,w,h;
207416022c9SBarry Smith   DrawCtx draw = (DrawCtx) viewer;
208416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
209416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
210416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
211416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
212416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
213416022c9SBarry Smith     yl = m - i - 1.0; yr = yl + 1.0;
214416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
215416022c9SBarry Smith       xl = a->j[j] + shift; xr = xl + 1.0;
216416022c9SBarry Smith       DrawRectangle(draw,xl,yl,xr,yr,DRAW_BLACK,DRAW_BLACK,DRAW_BLACK,DRAW_BLACK);
217416022c9SBarry Smith     }
218416022c9SBarry Smith   }
219416022c9SBarry Smith   DrawFlush(draw);
220416022c9SBarry Smith   return 0;
221416022c9SBarry Smith }
222416022c9SBarry Smith 
223416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
224416022c9SBarry Smith {
225416022c9SBarry Smith   Mat         A = (Mat) obj;
226416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
227416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
228416022c9SBarry Smith 
229416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
230416022c9SBarry Smith   if (!viewer) {
231416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
232416022c9SBarry Smith   }
233416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
234416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
235416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
236416022c9SBarry Smith     }
237416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
238416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
239416022c9SBarry Smith     }
240416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
241416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
242416022c9SBarry Smith     }
243416022c9SBarry Smith   }
244416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
245416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
246416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
24717ab2063SBarry Smith   }
24817ab2063SBarry Smith   return 0;
24917ab2063SBarry Smith }
25017ab2063SBarry Smith 
251416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
25217ab2063SBarry Smith {
253416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
254416022c9SBarry Smith   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
255416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
256416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
25717ab2063SBarry Smith 
25817ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
25917ab2063SBarry Smith 
26017ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
261416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
26217ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
26317ab2063SBarry Smith     if (fshift) {
264416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
26517ab2063SBarry Smith       N = ailen[i];
26617ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
26717ab2063SBarry Smith         ip[j-fshift] = ip[j];
26817ab2063SBarry Smith         ap[j-fshift] = ap[j];
26917ab2063SBarry Smith       }
27017ab2063SBarry Smith     }
27117ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
27217ab2063SBarry Smith   }
27317ab2063SBarry Smith   if (m) {
27417ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
27517ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
27617ab2063SBarry Smith   }
27717ab2063SBarry Smith   /* reset ilen and imax for each row */
27817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
27917ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
28017ab2063SBarry Smith   }
281416022c9SBarry Smith   a->nz = ai[m] + shift;
28217ab2063SBarry Smith 
28317ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
284416022c9SBarry Smith   if (fshift && a->diag) {
285416022c9SBarry Smith     PETSCFREE(a->diag);
286416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
287416022c9SBarry Smith     a->diag = 0;
28817ab2063SBarry Smith   }
289416022c9SBarry Smith   a->assembled = 1;
29017ab2063SBarry Smith   return 0;
29117ab2063SBarry Smith }
29217ab2063SBarry Smith 
293416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
29417ab2063SBarry Smith {
295416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
296416022c9SBarry Smith   PetscZero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
29717ab2063SBarry Smith   return 0;
29817ab2063SBarry Smith }
299416022c9SBarry Smith 
30017ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
30117ab2063SBarry Smith {
302416022c9SBarry Smith   Mat        A  = (Mat) obj;
303416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
30417ab2063SBarry Smith #if defined(PETSC_LOG)
305416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
30617ab2063SBarry Smith #endif
307416022c9SBarry Smith   PETSCFREE(a->a);
308416022c9SBarry Smith   if (!a->singlemalloc) { PETSCFREE(a->i); PETSCFREE(a->j);}
309416022c9SBarry Smith   if (a->diag) PETSCFREE(a->diag);
310416022c9SBarry Smith   if (a->ilen) PETSCFREE(a->ilen);
311416022c9SBarry Smith   if (a->imax) PETSCFREE(a->imax);
312416022c9SBarry Smith   if (a->solve_work) PETSCFREE(a->solve_work);
313416022c9SBarry Smith   PETSCFREE(a);
314416022c9SBarry Smith   PLogObjectDestroy(A);
315416022c9SBarry Smith   PETSCHEADERDESTROY(A);
31617ab2063SBarry Smith   return 0;
31717ab2063SBarry Smith }
31817ab2063SBarry Smith 
319416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
32017ab2063SBarry Smith {
32117ab2063SBarry Smith   return 0;
32217ab2063SBarry Smith }
32317ab2063SBarry Smith 
324416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
32517ab2063SBarry Smith {
326416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
327416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
328416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
329416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
330416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
331e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
332e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
333e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
334e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
335df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
336df719601SLois Curfman McInnes   else if (op == COLUMN_ORIENTED)
3374d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");}
338df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
3394d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
340e2f28af5SBarry Smith   else
3414d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
34217ab2063SBarry Smith   return 0;
34317ab2063SBarry Smith }
34417ab2063SBarry Smith 
345416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
34617ab2063SBarry Smith {
347416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
348416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
34917ab2063SBarry Smith   Scalar     *x, zero = 0.0;
35017ab2063SBarry Smith 
351416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
35217ab2063SBarry Smith   VecSet(&zero,v);
35317ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
354416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
355416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
356416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
357416022c9SBarry Smith       if (a->j[j]+shift == i) {
358416022c9SBarry Smith         x[i] = a->a[j];
35917ab2063SBarry Smith         break;
36017ab2063SBarry Smith       }
36117ab2063SBarry Smith     }
36217ab2063SBarry Smith   }
36317ab2063SBarry Smith   return 0;
36417ab2063SBarry Smith }
36517ab2063SBarry Smith 
36617ab2063SBarry Smith /* -------------------------------------------------------*/
36717ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
36817ab2063SBarry Smith /* -------------------------------------------------------*/
369416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
37017ab2063SBarry Smith {
371416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
37217ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
373416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
37417ab2063SBarry Smith 
375416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
37617ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
377416022c9SBarry Smith   PetscZero(y,a->n*sizeof(Scalar));
37817ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
37917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
380416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
381416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
382416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
38317ab2063SBarry Smith     alpha = x[i];
38417ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
38517ab2063SBarry Smith   }
386416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
38717ab2063SBarry Smith   return 0;
38817ab2063SBarry Smith }
389416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
39017ab2063SBarry Smith {
391416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
39217ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
393416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
39417ab2063SBarry Smith 
395416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
39617ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
39717ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
39817ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
39917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
400416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
401416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
402416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
40317ab2063SBarry Smith     alpha = x[i];
40417ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
40517ab2063SBarry Smith   }
40617ab2063SBarry Smith   return 0;
40717ab2063SBarry Smith }
40817ab2063SBarry Smith 
409416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
41017ab2063SBarry Smith {
411416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
41217ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
413416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
41417ab2063SBarry Smith 
415416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
41617ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
41717ab2063SBarry Smith   x = x + shift; /* shift for Fortran start by 1 indexing */
418416022c9SBarry Smith   idx  = a->j;
419416022c9SBarry Smith   v    = a->a;
420416022c9SBarry Smith   ii   = a->i;
42117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
422416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
42317ab2063SBarry Smith     sum  = 0.0;
42417ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
42517ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
426416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
42717ab2063SBarry Smith     y[i] = sum;
42817ab2063SBarry Smith   }
429416022c9SBarry Smith   PLogFlops(2*a->nz - m);
43017ab2063SBarry Smith   return 0;
43117ab2063SBarry Smith }
43217ab2063SBarry Smith 
433416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
43417ab2063SBarry Smith {
435416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
43617ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
437416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
43817ab2063SBarry Smith 
43948d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
44017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
44117ab2063SBarry Smith   x = x + shift; /* shift for Fortran start by 1 indexing */
44217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
443416022c9SBarry Smith     idx  = a->j + a->i[i] + shift;
444416022c9SBarry Smith     v    = a->a + a->i[i] + shift;
445416022c9SBarry Smith     n    = a->i[i+1] - a->i[i];
44617ab2063SBarry Smith     sum  = y[i];
44717ab2063SBarry Smith     SPARSEDENSEDOT(sum,x,v,idx,n);
44817ab2063SBarry Smith     z[i] = sum;
44917ab2063SBarry Smith   }
450416022c9SBarry Smith   PLogFlops(2*a->nz);
45117ab2063SBarry Smith   return 0;
45217ab2063SBarry Smith }
45317ab2063SBarry Smith 
45417ab2063SBarry Smith /*
45517ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
45617ab2063SBarry Smith */
45717ab2063SBarry Smith 
458416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
45917ab2063SBarry Smith {
460416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
461416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
46217ab2063SBarry Smith 
463416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
46417ab2063SBarry Smith   diag = (int *) PETSCMALLOC( (m+1)*sizeof(int)); CHKPTRQ(diag);
465416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
466416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
467416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
468416022c9SBarry Smith       if (a->j[j]+shift == i) {
46917ab2063SBarry Smith         diag[i] = j - shift;
47017ab2063SBarry Smith         break;
47117ab2063SBarry Smith       }
47217ab2063SBarry Smith     }
47317ab2063SBarry Smith   }
474416022c9SBarry Smith   a->diag = diag;
47517ab2063SBarry Smith   return 0;
47617ab2063SBarry Smith }
47717ab2063SBarry Smith 
478416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
47917ab2063SBarry Smith                            double fshift,int its,Vec xx)
48017ab2063SBarry Smith {
481416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
482416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
483416022c9SBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i;
484416022c9SBarry Smith   int        shift = a->indexshift;
48517ab2063SBarry Smith 
48617ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
487416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
488416022c9SBarry Smith   diag = a->diag;
489416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
49017ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
49117ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
49217ab2063SBarry Smith     bs = b + shift;
49317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
494416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
495416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
496416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
497416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
49817ab2063SBarry Smith         sum  = b[i]*d/omega;
49917ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
50017ab2063SBarry Smith         x[i] = sum;
50117ab2063SBarry Smith     }
50217ab2063SBarry Smith     return 0;
50317ab2063SBarry Smith   }
50417ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
50517ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
50617ab2063SBarry Smith   }
507416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
50817ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
50917ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
51017ab2063SBarry Smith 
51117ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
51217ab2063SBarry Smith 
51317ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
51417ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
51517ab2063SBarry Smith     is the relaxation factor.
51617ab2063SBarry Smith     */
51717ab2063SBarry Smith     t = (Scalar *) PETSCMALLOC( m*sizeof(Scalar) ); CHKPTRQ(t);
51817ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
51917ab2063SBarry Smith 
52017ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
52117ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
522416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
523416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
524416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
525416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
52617ab2063SBarry Smith       sum  = b[i];
52717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
52817ab2063SBarry Smith       x[i] = omega*(sum/d);
52917ab2063SBarry Smith     }
53017ab2063SBarry Smith 
53117ab2063SBarry Smith     /*  t = b - (2*E - D)x */
532416022c9SBarry Smith     v = a->a;
53317ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
53417ab2063SBarry Smith 
53517ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
536416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
537416022c9SBarry Smith     diag = a->diag;
53817ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
539416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
540416022c9SBarry Smith       n    = diag[i] - a->i[i];
541416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
542416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
54317ab2063SBarry Smith       sum  = t[i];
54417ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
54517ab2063SBarry Smith       t[i] = omega*(sum/d);
54617ab2063SBarry Smith     }
54717ab2063SBarry Smith 
54817ab2063SBarry Smith     /*  x = x + t */
54917ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
55017ab2063SBarry Smith     PETSCFREE(t);
55117ab2063SBarry Smith     return 0;
55217ab2063SBarry Smith   }
55317ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
55417ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
55517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
556416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
557416022c9SBarry Smith         n    = diag[i] - a->i[i];
558416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
559416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
56017ab2063SBarry Smith         sum  = b[i];
56117ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
56217ab2063SBarry Smith         x[i] = omega*(sum/d);
56317ab2063SBarry Smith       }
56417ab2063SBarry Smith       xb = x;
56517ab2063SBarry Smith     }
56617ab2063SBarry Smith     else xb = b;
56717ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
56817ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
56917ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
570416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
57117ab2063SBarry Smith       }
57217ab2063SBarry Smith     }
57317ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
57417ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
575416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
576416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
577416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
578416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
57917ab2063SBarry Smith         sum  = xb[i];
58017ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
58117ab2063SBarry Smith         x[i] = omega*(sum/d);
58217ab2063SBarry Smith       }
58317ab2063SBarry Smith     }
58417ab2063SBarry Smith     its--;
58517ab2063SBarry Smith   }
58617ab2063SBarry Smith   while (its--) {
58717ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
58817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
589416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
590416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
591416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
592416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
59317ab2063SBarry Smith         sum  = b[i];
59417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
59517ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
59617ab2063SBarry Smith       }
59717ab2063SBarry Smith     }
59817ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
59917ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
600416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
601416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
602416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
603416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
60417ab2063SBarry Smith         sum  = b[i];
60517ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
60617ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
60717ab2063SBarry Smith       }
60817ab2063SBarry Smith     }
60917ab2063SBarry Smith   }
61017ab2063SBarry Smith   return 0;
61117ab2063SBarry Smith }
61217ab2063SBarry Smith 
613416022c9SBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,
61417ab2063SBarry Smith                              int *nzalloc,int *mem)
61517ab2063SBarry Smith {
616416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
617416022c9SBarry Smith   *nz      = a->nz;
618416022c9SBarry Smith   *nzalloc = a->maxnz;
619416022c9SBarry Smith   *mem     = (int)A->mem;
62017ab2063SBarry Smith   return 0;
62117ab2063SBarry Smith }
62217ab2063SBarry Smith 
62317ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
62417ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
62517ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
62617ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
62717ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
62817ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
62917ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
63017ab2063SBarry Smith 
63117ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
63217ab2063SBarry Smith {
633416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
634416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
63517ab2063SBarry Smith 
63617ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
63717ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
63817ab2063SBarry Smith   if (diag) {
63917ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
640416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
641416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
642416022c9SBarry Smith         a->ilen[rows[i]] = 1;
643416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
644416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
64517ab2063SBarry Smith       }
64617ab2063SBarry Smith       else {
64717ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
64817ab2063SBarry Smith         CHKERRQ(ierr);
64917ab2063SBarry Smith       }
65017ab2063SBarry Smith     }
65117ab2063SBarry Smith   }
65217ab2063SBarry Smith   else {
65317ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
654416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
655416022c9SBarry Smith       a->ilen[rows[i]] = 0;
65617ab2063SBarry Smith     }
65717ab2063SBarry Smith   }
65817ab2063SBarry Smith   ISRestoreIndices(is,&rows);
65917ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
66017ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
66117ab2063SBarry Smith   return 0;
66217ab2063SBarry Smith }
66317ab2063SBarry Smith 
664416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
66517ab2063SBarry Smith {
666416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
667416022c9SBarry Smith   *m = a->m; *n = a->n;
66817ab2063SBarry Smith   return 0;
66917ab2063SBarry Smith }
67017ab2063SBarry Smith 
671416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
67217ab2063SBarry Smith {
673416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
674416022c9SBarry Smith   *m = 0; *n = a->m;
67517ab2063SBarry Smith   return 0;
67617ab2063SBarry Smith }
677416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
67817ab2063SBarry Smith {
679416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
680416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
68117ab2063SBarry Smith 
682416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
68317ab2063SBarry Smith 
684416022c9SBarry Smith   if (!a->assembled) {
685416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
686416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
68717ab2063SBarry Smith   }
688416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
689416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
69017ab2063SBarry Smith   if (idx) {
69117ab2063SBarry Smith     if (*nz) {
692416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
69317ab2063SBarry Smith       *idx = (int *) PETSCMALLOC( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
69417ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
69517ab2063SBarry Smith     }
69617ab2063SBarry Smith     else *idx = 0;
69717ab2063SBarry Smith   }
69817ab2063SBarry Smith   return 0;
69917ab2063SBarry Smith }
70017ab2063SBarry Smith 
701416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
70217ab2063SBarry Smith {
70317ab2063SBarry Smith   if (idx) {if (*idx) PETSCFREE(*idx);}
70417ab2063SBarry Smith   return 0;
70517ab2063SBarry Smith }
70617ab2063SBarry Smith 
707416022c9SBarry Smith static int MatNorm_SeqAIJ(Mat A,MatNormType type,double *norm)
70817ab2063SBarry Smith {
709416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
710416022c9SBarry Smith   Scalar     *v = a->a;
71117ab2063SBarry Smith   double     sum = 0.0;
712416022c9SBarry Smith   int        i, j,shift = a->indexshift;
71317ab2063SBarry Smith 
714416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
71517ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
716416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
71717ab2063SBarry Smith #if defined(PETSC_COMPLEX)
71817ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
71917ab2063SBarry Smith #else
72017ab2063SBarry Smith       sum += (*v)*(*v); v++;
72117ab2063SBarry Smith #endif
72217ab2063SBarry Smith     }
72317ab2063SBarry Smith     *norm = sqrt(sum);
72417ab2063SBarry Smith   }
72517ab2063SBarry Smith   else if (type == NORM_1) {
72617ab2063SBarry Smith     double *tmp;
727416022c9SBarry Smith     int    *jj = a->j;
728416022c9SBarry Smith     tmp = (double *) PETSCMALLOC( a->n*sizeof(double) ); CHKPTRQ(tmp);
729416022c9SBarry Smith     PetscZero(tmp,a->n*sizeof(double));
73017ab2063SBarry Smith     *norm = 0.0;
731416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
73217ab2063SBarry Smith #if defined(PETSC_COMPLEX)
73317ab2063SBarry Smith         tmp[*jj++ + shift] += abs(*v++);
73417ab2063SBarry Smith #else
73517ab2063SBarry Smith         tmp[*jj++ + shift] += fabs(*v++);
73617ab2063SBarry Smith #endif
73717ab2063SBarry Smith     }
738416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
73917ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
74017ab2063SBarry Smith     }
74117ab2063SBarry Smith     PETSCFREE(tmp);
74217ab2063SBarry Smith   }
74317ab2063SBarry Smith   else if (type == NORM_INFINITY) {
74417ab2063SBarry Smith     *norm = 0.0;
745416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
746416022c9SBarry Smith       v = a->a + a->i[j] + shift;
74717ab2063SBarry Smith       sum = 0.0;
748416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
74917ab2063SBarry Smith #if defined(PETSC_COMPLEX)
75017ab2063SBarry Smith         sum += abs(*v); v++;
75117ab2063SBarry Smith #else
75217ab2063SBarry Smith         sum += fabs(*v); v++;
75317ab2063SBarry Smith #endif
75417ab2063SBarry Smith       }
75517ab2063SBarry Smith       if (sum > *norm) *norm = sum;
75617ab2063SBarry Smith     }
75717ab2063SBarry Smith   }
75817ab2063SBarry Smith   else {
75948d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
76017ab2063SBarry Smith   }
76117ab2063SBarry Smith   return 0;
76217ab2063SBarry Smith }
76317ab2063SBarry Smith 
764416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
76517ab2063SBarry Smith {
766416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
767416022c9SBarry Smith   Mat        C;
768416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
769416022c9SBarry Smith   Scalar     *array = a->a;
770416022c9SBarry Smith   int        shift = a->indexshift;
77117ab2063SBarry Smith 
772416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
773416022c9SBarry Smith   col = (int *) PETSCMALLOC((1+a->n)*sizeof(int)); CHKPTRQ(col);
774416022c9SBarry Smith   PetscZero(col,(1+a->n)*sizeof(int));
77517ab2063SBarry Smith   if (shift) {
77617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
77717ab2063SBarry Smith   }
77817ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
779416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
78017ab2063SBarry Smith   PETSCFREE(col);
78117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
78217ab2063SBarry Smith     len = ai[i+1]-ai[i];
783416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
78417ab2063SBarry Smith     array += len; aj += len;
78517ab2063SBarry Smith   }
78617ab2063SBarry Smith   if (shift) {
78717ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
78817ab2063SBarry Smith   }
78917ab2063SBarry Smith 
790416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
791416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
79217ab2063SBarry Smith 
793416022c9SBarry Smith   if (B) {
794416022c9SBarry Smith     *B = C;
79517ab2063SBarry Smith   } else {
796416022c9SBarry Smith     /* This isn't really an in-place transpose */
797416022c9SBarry Smith     PETSCFREE(a->a);
798416022c9SBarry Smith     if (!a->singlemalloc) {PETSCFREE(a->i); PETSCFREE(a->j);}
799416022c9SBarry Smith     if (a->diag) PETSCFREE(a->diag);
800416022c9SBarry Smith     if (a->ilen) PETSCFREE(a->ilen);
801416022c9SBarry Smith     if (a->imax) PETSCFREE(a->imax);
802416022c9SBarry Smith     if (a->solve_work) PETSCFREE(a->solve_work);
803416022c9SBarry Smith     PETSCFREE(a);
804416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
805416022c9SBarry Smith     PETSCHEADERDESTROY(C);
80617ab2063SBarry Smith   }
80717ab2063SBarry Smith   return 0;
80817ab2063SBarry Smith }
80917ab2063SBarry Smith 
810416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
81117ab2063SBarry Smith {
812416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
81317ab2063SBarry Smith   Scalar     *l,*r,x,*v;
814416022c9SBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj;
815416022c9SBarry Smith   int        shift = a->indexshift;
81617ab2063SBarry Smith 
81748d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
81817ab2063SBarry Smith   if (ll) {
81917ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
820416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
821416022c9SBarry Smith     v = a->a;
82217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
82317ab2063SBarry Smith       x = l[i];
824416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
82517ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
82617ab2063SBarry Smith     }
82717ab2063SBarry Smith   }
82817ab2063SBarry Smith   if (rr) {
82917ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
830416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
831416022c9SBarry Smith     v = a->a; jj = a->j;
83217ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
83317ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
83417ab2063SBarry Smith     }
83517ab2063SBarry Smith   }
83617ab2063SBarry Smith   return 0;
83717ab2063SBarry Smith }
83817ab2063SBarry Smith 
839416022c9SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,Mat *B)
84017ab2063SBarry Smith {
841db02288aSLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c;
84202834360SBarry Smith   int        nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
84302834360SBarry Smith   int        *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
844db02288aSLois Curfman McInnes   int        first,step,*starts,*j_new,*i_new;
845db02288aSLois Curfman McInnes   Scalar     *vwork,*a_new;
846416022c9SBarry Smith   Mat        C;
84717ab2063SBarry Smith 
848416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
84917ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
85017ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
85117ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
85217ab2063SBarry Smith 
85302834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
85402834360SBarry Smith     /* special case of contiguous rows */
855db02288aSLois Curfman McInnes     lens   = (int *) PETSCMALLOC((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
85602834360SBarry Smith     starts = lens + ncols;
85702834360SBarry Smith     /* loop over new rows determining lens and starting points */
85802834360SBarry Smith     for (i=0; i<nrows; i++) {
85902834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
86002834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
86102834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
86202834360SBarry Smith         if (a->j[k] >= first) {
86302834360SBarry Smith           starts[i] = k;
86402834360SBarry Smith           break;
86502834360SBarry Smith 	}
86602834360SBarry Smith       }
86702834360SBarry Smith       lens[i] = 0;
86802834360SBarry Smith       while (k < kend) {
86902834360SBarry Smith         if (a->j[k++] >= first+ncols) break;
87002834360SBarry Smith         lens[i]++;
87102834360SBarry Smith       }
87202834360SBarry Smith     }
87302834360SBarry Smith     /* create submatrix */
874*08480c60SBarry Smith     if (MatValidMatrix(*B)) {
875*08480c60SBarry Smith       int n_cols,n_rows;
876*08480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
877*08480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
878*08480c60SBarry Smith       MatZeroEntries(*B);
879*08480c60SBarry Smith       C = *B;
880*08480c60SBarry Smith     }
881*08480c60SBarry Smith     else {
88202834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
883*08480c60SBarry Smith     }
884db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
885db02288aSLois Curfman McInnes 
88602834360SBarry Smith     /* loop over rows inserting into submatrix */
887db02288aSLois Curfman McInnes     a_new    = c->a;
888db02288aSLois Curfman McInnes     j_new    = c->j;
889db02288aSLois Curfman McInnes     i_new    = c->i;
890db02288aSLois Curfman McInnes     i_new[0] = -shift;
89102834360SBarry Smith     for (i=0; i<nrows; i++) {
89202834360SBarry Smith       for ( k=0; k<lens[i]; k++ ) {
893db02288aSLois Curfman McInnes         *j_new++ = a->j[k+starts[i]] - first;
89402834360SBarry Smith       }
895db02288aSLois Curfman McInnes       PetscMemcpy(a_new,a->a + starts[i],lens[i]*sizeof(Scalar));
896db02288aSLois Curfman McInnes       a_new      += lens[i];
897db02288aSLois Curfman McInnes       i_new[i+1]  = i_new[i] + lens[i];
8981987afe7SBarry Smith       c->ilen[i]  = lens[i];
89902834360SBarry Smith     }
900db02288aSLois Curfman McInnes     PETSCFREE(lens);
90102834360SBarry Smith   }
90202834360SBarry Smith   else {
90302834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
90417ab2063SBarry Smith     smap  = (int *) PETSCMALLOC((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
90502834360SBarry Smith     ssmap = smap + shift;
90602834360SBarry Smith     cwork = (int *) PETSCMALLOC((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork);
90702834360SBarry Smith     lens  = cwork + ncols;
90817ab2063SBarry Smith     vwork = (Scalar *) PETSCMALLOC((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
909416022c9SBarry Smith     PetscZero(smap,oldcols*sizeof(int));
91017ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
91102834360SBarry Smith     /* determine lens of each row */
91202834360SBarry Smith     for (i=0; i<nrows; i++) {
91302834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
91402834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
91502834360SBarry Smith       lens[i] = 0;
91602834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
91702834360SBarry Smith         if (ssmap[a->j[k]]) {
91802834360SBarry Smith           lens[i]++;
91902834360SBarry Smith         }
92002834360SBarry Smith       }
92102834360SBarry Smith     }
92217ab2063SBarry Smith     /* Create and fill new matrix */
923*08480c60SBarry Smith     if (MatValidMatrix(*B)) {
924*08480c60SBarry Smith       int n_cols,n_rows;
925*08480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
926*08480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
927*08480c60SBarry Smith       MatZeroEntries(*B);
928*08480c60SBarry Smith       C = *B;
929*08480c60SBarry Smith     }
930*08480c60SBarry Smith     else {
93102834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
932*08480c60SBarry Smith     }
93317ab2063SBarry Smith     for (i=0; i<nrows; i++) {
93417ab2063SBarry Smith       nznew  = 0;
935416022c9SBarry Smith       kstart = a->i[irow[i]]+shift;
936416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
93717ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
93802834360SBarry Smith         if (ssmap[a->j[k]]) {
93902834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
940416022c9SBarry Smith           vwork[nznew++] = a->a[k];
94117ab2063SBarry Smith         }
94217ab2063SBarry Smith       }
943416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
94417ab2063SBarry Smith     }
94502834360SBarry Smith     /* Free work space */
94602834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
94702834360SBarry Smith     PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
94802834360SBarry Smith   }
949416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
950416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
95117ab2063SBarry Smith 
95217ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
953416022c9SBarry Smith   *B = C;
95417ab2063SBarry Smith   return 0;
95517ab2063SBarry Smith }
95617ab2063SBarry Smith 
957a871dcd8SBarry Smith /*
95863b91edcSBarry Smith      note: This can only work for identity for row and col. It would
95963b91edcSBarry Smith    be good to check this and otherwise generate an error.
960a871dcd8SBarry Smith */
96163b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
962a871dcd8SBarry Smith {
96363b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
964*08480c60SBarry Smith   int        ierr;
96563b91edcSBarry Smith   Mat        outA;
96663b91edcSBarry Smith 
967a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
968a871dcd8SBarry Smith 
96963b91edcSBarry Smith   outA          = inA;
97063b91edcSBarry Smith   inA->factor   = FACTOR_LU;
97163b91edcSBarry Smith   a->row        = row;
97263b91edcSBarry Smith   a->col        = col;
97363b91edcSBarry Smith 
97463b91edcSBarry Smith   a->solve_work = (Scalar *) PETSCMALLOC( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
97563b91edcSBarry Smith 
976*08480c60SBarry Smith   if (!a->diag) {
977*08480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
97863b91edcSBarry Smith   }
97963b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
980a871dcd8SBarry Smith   return 0;
981a871dcd8SBarry Smith }
982a871dcd8SBarry Smith 
98317ab2063SBarry Smith /* -------------------------------------------------------------------*/
98417ab2063SBarry Smith 
98517ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
98617ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
987416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
988416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
98917ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
99017ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
99117ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
99217ab2063SBarry Smith        MatRelax_SeqAIJ,
99317ab2063SBarry Smith        MatTranspose_SeqAIJ,
99417ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
99517ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
99617ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
99717ab2063SBarry Smith        MatCompress_SeqAIJ,
99817ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
99917ab2063SBarry Smith        MatGetReordering_SeqAIJ,
100017ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
100117ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
100217ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
100317ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1004416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
1005a871dcd8SBarry Smith        MatCopyPrivate_SeqAIJ,0,0,
100663b91edcSBarry Smith        MatILUFactor_SeqAIJ};
100717ab2063SBarry Smith 
100817ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
100917ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
101017ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
101117ab2063SBarry Smith 
101217ab2063SBarry Smith /*@C
101317ab2063SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
101417ab2063SBarry Smith    (the default uniprocessor PETSc format).
101517ab2063SBarry Smith 
101617ab2063SBarry Smith    Input Parameters:
101717ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
101817ab2063SBarry Smith .  m - number of rows
101917ab2063SBarry Smith .  n - number of columns
102017ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
102117ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
102217ab2063SBarry Smith 
102317ab2063SBarry Smith    Output Parameter:
1024416022c9SBarry Smith .  A - the matrix
102517ab2063SBarry Smith 
102617ab2063SBarry Smith    Notes:
102717ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
102817ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
10290002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
10300002213bSLois Curfman McInnes    either one (as in Fortran) or zero.  See the users manual for details.
103117ab2063SBarry Smith 
103217ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
103317ab2063SBarry Smith    Set both nz and nnz to zero for PETSc to control dynamic memory
103417ab2063SBarry Smith    allocation.
103517ab2063SBarry Smith 
103617ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse
103717ab2063SBarry Smith 
103817ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
103917ab2063SBarry Smith @*/
1040416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
104117ab2063SBarry Smith {
1042416022c9SBarry Smith   Mat        B;
1043416022c9SBarry Smith   Mat_SeqAIJ *b;
104417ab2063SBarry Smith   int        i,len,ierr;
1045416022c9SBarry Smith   *A      = 0;
1046416022c9SBarry Smith   PETSCHEADERCREATE(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1047416022c9SBarry Smith   PLogObjectCreate(B);
1048416022c9SBarry Smith   B->data               = (void *) (b = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(b);
1049416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1050416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1051416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1052416022c9SBarry Smith   B->factor           = 0;
1053416022c9SBarry Smith   B->lupivotthreshold = 1.0;
1054416022c9SBarry Smith   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1055416022c9SBarry Smith   b->row              = 0;
1056416022c9SBarry Smith   b->col              = 0;
1057416022c9SBarry Smith   b->indexshift       = 0;
1058416022c9SBarry Smith   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
105917ab2063SBarry Smith 
1060416022c9SBarry Smith   b->m       = m;
1061416022c9SBarry Smith   b->n       = n;
1062416022c9SBarry Smith   b->imax    = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
106317ab2063SBarry Smith   if (!nnz) {
106417ab2063SBarry Smith     if (nz <= 0) nz = 1;
1065416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
106617ab2063SBarry Smith     nz = nz*m;
106717ab2063SBarry Smith   }
106817ab2063SBarry Smith   else {
106917ab2063SBarry Smith     nz = 0;
1070416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
107117ab2063SBarry Smith   }
107217ab2063SBarry Smith 
107317ab2063SBarry Smith   /* allocate the matrix space */
1074416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
1075416022c9SBarry Smith   b->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(b->a);
1076416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1077416022c9SBarry Smith   PetscZero(b->j,nz*sizeof(int));
1078416022c9SBarry Smith   b->i  = b->j + nz;
1079416022c9SBarry Smith   b->singlemalloc = 1;
108017ab2063SBarry Smith 
1081416022c9SBarry Smith   b->i[0] = -b->indexshift;
108217ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1083416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
108417ab2063SBarry Smith   }
108517ab2063SBarry Smith 
1086416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
1087416022c9SBarry Smith   b->ilen = (int *) PETSCMALLOC((m+1)*sizeof(int));
1088416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1089416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
109017ab2063SBarry Smith 
1091416022c9SBarry Smith   b->nz          = 0;
1092416022c9SBarry Smith   b->maxnz       = nz;
1093416022c9SBarry Smith   b->sorted      = 0;
1094416022c9SBarry Smith   b->roworiented = 1;
1095416022c9SBarry Smith   b->nonew       = 0;
1096416022c9SBarry Smith   b->diag        = 0;
1097416022c9SBarry Smith   b->assembled   = 0;
1098416022c9SBarry Smith   b->solve_work  = 0;
109971bd300dSLois Curfman McInnes   b->spptr       = 0;
110017ab2063SBarry Smith 
1101416022c9SBarry Smith   *A = B;
110217ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_superlu")) {
1103416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
110417ab2063SBarry Smith   }
110517ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_essl")) {
1106416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
110717ab2063SBarry Smith   }
110817ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_dxml")) {
1109416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1110416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
111117ab2063SBarry Smith   }
111217ab2063SBarry Smith 
111317ab2063SBarry Smith   return 0;
111417ab2063SBarry Smith }
111517ab2063SBarry Smith 
1116*08480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues)
111717ab2063SBarry Smith {
1118416022c9SBarry Smith   Mat        C;
1119416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1120*08480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
112117ab2063SBarry Smith 
11224043dd9cSLois Curfman McInnes   *B = 0;
1123416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
1124416022c9SBarry Smith   PETSCHEADERCREATE(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1125416022c9SBarry Smith   PLogObjectCreate(C);
1126416022c9SBarry Smith   C->data       = (void *) (c = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(c);
1127416022c9SBarry Smith   PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps));
1128416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1129416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1130416022c9SBarry Smith   C->factor     = A->factor;
1131416022c9SBarry Smith   c->row        = 0;
1132416022c9SBarry Smith   c->col        = 0;
1133416022c9SBarry Smith   c->indexshift = shift;
113417ab2063SBarry Smith 
1135416022c9SBarry Smith   c->m          = a->m;
1136416022c9SBarry Smith   c->n          = a->n;
113717ab2063SBarry Smith 
1138416022c9SBarry Smith   c->imax       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1139416022c9SBarry Smith   c->ilen       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
114017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1141416022c9SBarry Smith     c->imax[i] = a->imax[i];
1142416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
114317ab2063SBarry Smith   }
114417ab2063SBarry Smith 
114517ab2063SBarry Smith   /* allocate the matrix space */
1146416022c9SBarry Smith   c->singlemalloc = 1;
1147416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1148416022c9SBarry Smith   c->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(c->a);
1149416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1150416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1151416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
115217ab2063SBarry Smith   if (m > 0) {
1153416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1154*08480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1155416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
115617ab2063SBarry Smith     }
1157*08480c60SBarry Smith   }
115817ab2063SBarry Smith 
1159416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1160416022c9SBarry Smith   c->sorted      = a->sorted;
1161416022c9SBarry Smith   c->roworiented = a->roworiented;
1162416022c9SBarry Smith   c->nonew       = a->nonew;
116317ab2063SBarry Smith 
1164416022c9SBarry Smith   if (a->diag) {
1165416022c9SBarry Smith     c->diag = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1166416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
116717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1168416022c9SBarry Smith       c->diag[i] = a->diag[i];
116917ab2063SBarry Smith     }
117017ab2063SBarry Smith   }
1171416022c9SBarry Smith   else c->diag        = 0;
1172416022c9SBarry Smith   c->assembled        = 1;
1173416022c9SBarry Smith   c->nz               = a->nz;
1174416022c9SBarry Smith   c->maxnz            = a->maxnz;
1175416022c9SBarry Smith   c->solve_work       = 0;
11764043dd9cSLois Curfman McInnes   c->spptr            = 0;
1177416022c9SBarry Smith   *B = C;
117817ab2063SBarry Smith   return 0;
117917ab2063SBarry Smith }
118017ab2063SBarry Smith 
1181416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
118217ab2063SBarry Smith {
1183416022c9SBarry Smith   Mat_SeqAIJ   *a;
1184416022c9SBarry Smith   Mat          B;
118517699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
118617ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
118717ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
118817ab2063SBarry Smith 
118917699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
119017699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
119117ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1192416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
119348d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
119417ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
119517ab2063SBarry Smith 
119617ab2063SBarry Smith   /* read in row lengths */
119717ab2063SBarry Smith   rowlengths = (int*) PETSCMALLOC( M*sizeof(int) ); CHKPTRQ(rowlengths);
1198416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
119917ab2063SBarry Smith 
120017ab2063SBarry Smith   /* create our matrix */
1201416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1202416022c9SBarry Smith   B = *A;
1203416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1204416022c9SBarry Smith   shift = a->indexshift;
120517ab2063SBarry Smith 
120617ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1207416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
120817ab2063SBarry Smith   if (shift) {
120917ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1210416022c9SBarry Smith       a->j[i] += 1;
121117ab2063SBarry Smith     }
121217ab2063SBarry Smith   }
121317ab2063SBarry Smith 
121417ab2063SBarry Smith   /* read in nonzero values */
1215416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
121617ab2063SBarry Smith 
121717ab2063SBarry Smith   /* set matrix "i" values */
1218416022c9SBarry Smith   a->i[0] = -shift;
121917ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1220416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1221416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
122217ab2063SBarry Smith   }
122317ab2063SBarry Smith   PETSCFREE(rowlengths);
122417ab2063SBarry Smith 
1225416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1226416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
122717ab2063SBarry Smith   return 0;
122817ab2063SBarry Smith }
122917ab2063SBarry Smith 
123017ab2063SBarry Smith 
123117ab2063SBarry Smith 
1232