xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 48d914877651b9fa76f196ebaf17c108719d05c1)
117ab2063SBarry Smith #ifndef lint
2*48d91487SBarry Smith static char vcid[] = "$Id: aij.c,v 1.92 1995/09/30 19:28:44 bsmith 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);
120416022c9SBarry Smith 
121416022c9SBarry Smith   col_lens    = (int *) PETSCMALLOC( (4+a->nz)*sizeof(int) ); CHKPTRQ(col_lens);
122416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
123416022c9SBarry Smith   col_lens[1] = a->m;
124416022c9SBarry Smith   col_lens[2] = a->n;
125416022c9SBarry Smith   col_lens[3] = a->nz;
126416022c9SBarry Smith 
127416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
128416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
129416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
13017ab2063SBarry Smith   }
131416022c9SBarry Smith   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
132416022c9SBarry Smith   PETSCFREE(col_lens);
133416022c9SBarry Smith 
134416022c9SBarry Smith   /* store column indices (zero start index) */
135416022c9SBarry Smith   if (a->indexshift) {
136416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
13717ab2063SBarry Smith   }
138416022c9SBarry Smith   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
139416022c9SBarry Smith   if (a->indexshift) {
140416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
14117ab2063SBarry Smith   }
142416022c9SBarry Smith 
143416022c9SBarry Smith   /* store nonzero values */
144416022c9SBarry Smith   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
14517ab2063SBarry Smith   return 0;
14617ab2063SBarry Smith }
147416022c9SBarry Smith 
148416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
149416022c9SBarry Smith {
150416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
151416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
15217ab2063SBarry Smith   FILE        *fd;
15317ab2063SBarry Smith   char        *outputname;
15417ab2063SBarry Smith 
155416022c9SBarry Smith   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
156416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
157416022c9SBarry Smith   ierr = ViewerFileGetFormat_Private(viewer,&format);
15817ab2063SBarry Smith   if (format == FILE_FORMAT_INFO) {
159416022c9SBarry Smith     ; /* do nothing for now */
16017ab2063SBarry Smith   }
16117ab2063SBarry Smith   else if (format == FILE_FORMAT_MATLAB) {
16217ab2063SBarry Smith     int nz, nzalloc, mem;
163416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
164416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
16517ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
16617ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
16717ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
16817ab2063SBarry Smith 
16917ab2063SBarry Smith     for (i=0; i<m; i++) {
170416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
17117ab2063SBarry Smith #if defined(PETSC_COMPLEX)
172416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
173416022c9SBarry Smith                    imag(a->a[j]));
17417ab2063SBarry Smith #else
175416022c9SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
17617ab2063SBarry Smith #endif
17717ab2063SBarry Smith       }
17817ab2063SBarry Smith     }
17917ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
18017ab2063SBarry Smith   }
18117ab2063SBarry Smith   else {
18217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
18317ab2063SBarry Smith       fprintf(fd,"row %d:",i);
184416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
18517ab2063SBarry Smith #if defined(PETSC_COMPLEX)
186416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
187416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
18817ab2063SBarry Smith         }
18917ab2063SBarry Smith         else {
190416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
19117ab2063SBarry Smith         }
19217ab2063SBarry Smith #else
193416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
19417ab2063SBarry Smith #endif
19517ab2063SBarry Smith       }
19617ab2063SBarry Smith       fprintf(fd,"\n");
19717ab2063SBarry Smith     }
19817ab2063SBarry Smith   }
19917ab2063SBarry Smith   fflush(fd);
200416022c9SBarry Smith   return 0;
201416022c9SBarry Smith }
202416022c9SBarry Smith 
203416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
204416022c9SBarry Smith {
205416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
206416022c9SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift;
207416022c9SBarry Smith   double      xl,yl,xr,yr,w,h;
208416022c9SBarry Smith   DrawCtx draw = (DrawCtx) viewer;
209416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
210416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
211416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
212416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
213416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
214416022c9SBarry Smith     yl = m - i - 1.0; yr = yl + 1.0;
215416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
216416022c9SBarry Smith       xl = a->j[j] + shift; xr = xl + 1.0;
217416022c9SBarry Smith       DrawRectangle(draw,xl,yl,xr,yr,DRAW_BLACK,DRAW_BLACK,DRAW_BLACK,DRAW_BLACK);
218416022c9SBarry Smith     }
219416022c9SBarry Smith   }
220416022c9SBarry Smith   DrawFlush(draw);
221416022c9SBarry Smith   return 0;
222416022c9SBarry Smith }
223416022c9SBarry Smith 
224416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
225416022c9SBarry Smith {
226416022c9SBarry Smith   Mat         A = (Mat) obj;
227416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
228416022c9SBarry Smith   PetscObject vobj = (PetscObject) viewer;
229416022c9SBarry Smith 
230416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
231416022c9SBarry Smith   if (!viewer) {
232416022c9SBarry Smith     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
233416022c9SBarry Smith   }
234416022c9SBarry Smith   if (vobj->cookie == VIEWER_COOKIE) {
235416022c9SBarry Smith     if (vobj->type == MATLAB_VIEWER) {
236416022c9SBarry Smith       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
237416022c9SBarry Smith     }
238416022c9SBarry Smith     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
239416022c9SBarry Smith       return MatView_SeqAIJ_ASCII(A,viewer);
240416022c9SBarry Smith     }
241416022c9SBarry Smith     else if (vobj->type == BINARY_FILE_VIEWER) {
242416022c9SBarry Smith       return MatView_SeqAIJ_Binary(A,viewer);
243416022c9SBarry Smith     }
244416022c9SBarry Smith   }
245416022c9SBarry Smith   else if (vobj->cookie == DRAW_COOKIE) {
246416022c9SBarry Smith     if (vobj->type == NULLWINDOW) return 0;
247416022c9SBarry Smith     else return MatView_SeqAIJ_Draw(A,viewer);
24817ab2063SBarry Smith   }
24917ab2063SBarry Smith   return 0;
25017ab2063SBarry Smith }
25117ab2063SBarry Smith 
252416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
25317ab2063SBarry Smith {
254416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
255416022c9SBarry Smith   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
256416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
257416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
25817ab2063SBarry Smith 
25917ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
26017ab2063SBarry Smith 
26117ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
262416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
26317ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
26417ab2063SBarry Smith     if (fshift) {
265416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
26617ab2063SBarry Smith       N = ailen[i];
26717ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
26817ab2063SBarry Smith         ip[j-fshift] = ip[j];
26917ab2063SBarry Smith         ap[j-fshift] = ap[j];
27017ab2063SBarry Smith       }
27117ab2063SBarry Smith     }
27217ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
27317ab2063SBarry Smith   }
27417ab2063SBarry Smith   if (m) {
27517ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
27617ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
27717ab2063SBarry Smith   }
27817ab2063SBarry Smith   /* reset ilen and imax for each row */
27917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
28017ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
28117ab2063SBarry Smith   }
282416022c9SBarry Smith   a->nz = ai[m] + shift;
28317ab2063SBarry Smith 
28417ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
285416022c9SBarry Smith   if (fshift && a->diag) {
286416022c9SBarry Smith     PETSCFREE(a->diag);
287416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
288416022c9SBarry Smith     a->diag = 0;
28917ab2063SBarry Smith   }
290416022c9SBarry Smith   a->assembled = 1;
29117ab2063SBarry Smith   return 0;
29217ab2063SBarry Smith }
29317ab2063SBarry Smith 
294416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
29517ab2063SBarry Smith {
296416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
297416022c9SBarry Smith   PetscZero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
29817ab2063SBarry Smith   return 0;
29917ab2063SBarry Smith }
300416022c9SBarry Smith 
30117ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
30217ab2063SBarry Smith {
303416022c9SBarry Smith   Mat        A  = (Mat) obj;
304416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
30517ab2063SBarry Smith #if defined(PETSC_LOG)
306416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
30717ab2063SBarry Smith #endif
308416022c9SBarry Smith   PETSCFREE(a->a);
309416022c9SBarry Smith   if (!a->singlemalloc) { PETSCFREE(a->i); PETSCFREE(a->j);}
310416022c9SBarry Smith   if (a->diag) PETSCFREE(a->diag);
311416022c9SBarry Smith   if (a->ilen) PETSCFREE(a->ilen);
312416022c9SBarry Smith   if (a->imax) PETSCFREE(a->imax);
313416022c9SBarry Smith   if (a->solve_work) PETSCFREE(a->solve_work);
314416022c9SBarry Smith   PETSCFREE(a);
315416022c9SBarry Smith   PLogObjectDestroy(A);
316416022c9SBarry Smith   PETSCHEADERDESTROY(A);
31717ab2063SBarry Smith   return 0;
31817ab2063SBarry Smith }
31917ab2063SBarry Smith 
320416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
32117ab2063SBarry Smith {
32217ab2063SBarry Smith   return 0;
32317ab2063SBarry Smith }
32417ab2063SBarry Smith 
325416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
32617ab2063SBarry Smith {
327416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
328416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
329416022c9SBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
330416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
33117ab2063SBarry Smith   /* doesn't care about sorted rows */
332416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
333416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
33417ab2063SBarry Smith 
335416022c9SBarry Smith   if (op == COLUMN_ORIENTED) SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:Column input not supported");
33617ab2063SBarry Smith   return 0;
33717ab2063SBarry Smith }
33817ab2063SBarry Smith 
339416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
34017ab2063SBarry Smith {
341416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
342416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
34317ab2063SBarry Smith   Scalar     *x, zero = 0.0;
34417ab2063SBarry Smith 
345416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
34617ab2063SBarry Smith   VecSet(&zero,v);
34717ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
348416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
349416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
350416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
351416022c9SBarry Smith       if (a->j[j]+shift == i) {
352416022c9SBarry Smith         x[i] = a->a[j];
35317ab2063SBarry Smith         break;
35417ab2063SBarry Smith       }
35517ab2063SBarry Smith     }
35617ab2063SBarry Smith   }
35717ab2063SBarry Smith   return 0;
35817ab2063SBarry Smith }
35917ab2063SBarry Smith 
36017ab2063SBarry Smith /* -------------------------------------------------------*/
36117ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
36217ab2063SBarry Smith /* -------------------------------------------------------*/
363416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
36417ab2063SBarry Smith {
365416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
36617ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
367416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
36817ab2063SBarry Smith 
369416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
37017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
371416022c9SBarry Smith   PetscZero(y,a->n*sizeof(Scalar));
37217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
37317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
374416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
375416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
376416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
37717ab2063SBarry Smith     alpha = x[i];
37817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
37917ab2063SBarry Smith   }
380416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
38117ab2063SBarry Smith   return 0;
38217ab2063SBarry Smith }
383416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
38417ab2063SBarry Smith {
385416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
38617ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
387416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
38817ab2063SBarry Smith 
389416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
39017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
39117ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
39217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
39317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
394416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
395416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
396416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
39717ab2063SBarry Smith     alpha = x[i];
39817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
39917ab2063SBarry Smith   }
40017ab2063SBarry Smith   return 0;
40117ab2063SBarry Smith }
40217ab2063SBarry Smith 
403416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
40417ab2063SBarry Smith {
405416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
40617ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
407416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
40817ab2063SBarry Smith 
409416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
41017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
41117ab2063SBarry Smith   x = x + shift; /* shift for Fortran start by 1 indexing */
412416022c9SBarry Smith   idx  = a->j;
413416022c9SBarry Smith   v    = a->a;
414416022c9SBarry Smith   ii   = a->i;
415416022c9SBarry Smith #if defined(PARCH_rs6000)
416416022c9SBarry Smith #pragma disjoint (*x,*v,*y)
417416022c9SBarry Smith #endif
41817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
419416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
42017ab2063SBarry Smith     sum  = 0.0;
42117ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
42217ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
423416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
42417ab2063SBarry Smith     y[i] = sum;
42517ab2063SBarry Smith   }
426416022c9SBarry Smith   PLogFlops(2*a->nz - m);
42717ab2063SBarry Smith   return 0;
42817ab2063SBarry Smith }
42917ab2063SBarry Smith 
430416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
43117ab2063SBarry Smith {
432416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
43317ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
434416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
43517ab2063SBarry Smith 
436*48d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
43717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
43817ab2063SBarry Smith   x = x + shift; /* shift for Fortran start by 1 indexing */
43917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
440416022c9SBarry Smith     idx  = a->j + a->i[i] + shift;
441416022c9SBarry Smith     v    = a->a + a->i[i] + shift;
442416022c9SBarry Smith     n    = a->i[i+1] - a->i[i];
44317ab2063SBarry Smith     sum  = y[i];
44417ab2063SBarry Smith     SPARSEDENSEDOT(sum,x,v,idx,n);
44517ab2063SBarry Smith     z[i] = sum;
44617ab2063SBarry Smith   }
447416022c9SBarry Smith   PLogFlops(2*a->nz);
44817ab2063SBarry Smith   return 0;
44917ab2063SBarry Smith }
45017ab2063SBarry Smith 
45117ab2063SBarry Smith /*
45217ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
45317ab2063SBarry Smith */
45417ab2063SBarry Smith 
455416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
45617ab2063SBarry Smith {
457416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
458416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
45917ab2063SBarry Smith 
460416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
46117ab2063SBarry Smith   diag = (int *) PETSCMALLOC( (m+1)*sizeof(int)); CHKPTRQ(diag);
462416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
463416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
464416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
465416022c9SBarry Smith       if (a->j[j]+shift == i) {
46617ab2063SBarry Smith         diag[i] = j - shift;
46717ab2063SBarry Smith         break;
46817ab2063SBarry Smith       }
46917ab2063SBarry Smith     }
47017ab2063SBarry Smith   }
471416022c9SBarry Smith   a->diag = diag;
47217ab2063SBarry Smith   return 0;
47317ab2063SBarry Smith }
47417ab2063SBarry Smith 
475416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
47617ab2063SBarry Smith                            double fshift,int its,Vec xx)
47717ab2063SBarry Smith {
478416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
479416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
480416022c9SBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i;
481416022c9SBarry Smith   int        shift = a->indexshift;
48217ab2063SBarry Smith 
48317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
484416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
485416022c9SBarry Smith   diag = a->diag;
486416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
48717ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
48817ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
48917ab2063SBarry Smith     bs = b + shift;
49017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
491416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
492416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
493416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
494416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
49517ab2063SBarry Smith         sum  = b[i]*d/omega;
49617ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
49717ab2063SBarry Smith         x[i] = sum;
49817ab2063SBarry Smith     }
49917ab2063SBarry Smith     return 0;
50017ab2063SBarry Smith   }
50117ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
50217ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
50317ab2063SBarry Smith   }
504416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
50517ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
50617ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
50717ab2063SBarry Smith 
50817ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
50917ab2063SBarry Smith 
51017ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
51117ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
51217ab2063SBarry Smith     is the relaxation factor.
51317ab2063SBarry Smith     */
51417ab2063SBarry Smith     t = (Scalar *) PETSCMALLOC( m*sizeof(Scalar) ); CHKPTRQ(t);
51517ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
51617ab2063SBarry Smith 
51717ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
51817ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
519416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
520416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
521416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
522416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
52317ab2063SBarry Smith       sum  = b[i];
52417ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
52517ab2063SBarry Smith       x[i] = omega*(sum/d);
52617ab2063SBarry Smith     }
52717ab2063SBarry Smith 
52817ab2063SBarry Smith     /*  t = b - (2*E - D)x */
529416022c9SBarry Smith     v = a->a;
53017ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
53117ab2063SBarry Smith 
53217ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
533416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
534416022c9SBarry Smith     diag = a->diag;
53517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
536416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
537416022c9SBarry Smith       n    = diag[i] - a->i[i];
538416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
539416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
54017ab2063SBarry Smith       sum  = t[i];
54117ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
54217ab2063SBarry Smith       t[i] = omega*(sum/d);
54317ab2063SBarry Smith     }
54417ab2063SBarry Smith 
54517ab2063SBarry Smith     /*  x = x + t */
54617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
54717ab2063SBarry Smith     PETSCFREE(t);
54817ab2063SBarry Smith     return 0;
54917ab2063SBarry Smith   }
55017ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
55117ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
55217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
553416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
554416022c9SBarry Smith         n    = diag[i] - a->i[i];
555416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
556416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
55717ab2063SBarry Smith         sum  = b[i];
55817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
55917ab2063SBarry Smith         x[i] = omega*(sum/d);
56017ab2063SBarry Smith       }
56117ab2063SBarry Smith       xb = x;
56217ab2063SBarry Smith     }
56317ab2063SBarry Smith     else xb = b;
56417ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
56517ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
56617ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
567416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
56817ab2063SBarry Smith       }
56917ab2063SBarry Smith     }
57017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
57117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
572416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
573416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
574416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
575416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
57617ab2063SBarry Smith         sum  = xb[i];
57717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
57817ab2063SBarry Smith         x[i] = omega*(sum/d);
57917ab2063SBarry Smith       }
58017ab2063SBarry Smith     }
58117ab2063SBarry Smith     its--;
58217ab2063SBarry Smith   }
58317ab2063SBarry Smith   while (its--) {
58417ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
58517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
586416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
587416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
588416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
589416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
59017ab2063SBarry Smith         sum  = b[i];
59117ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
59217ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
59317ab2063SBarry Smith       }
59417ab2063SBarry Smith     }
59517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
59617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
597416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
598416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
599416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
600416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
60117ab2063SBarry Smith         sum  = b[i];
60217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
60317ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
60417ab2063SBarry Smith       }
60517ab2063SBarry Smith     }
60617ab2063SBarry Smith   }
60717ab2063SBarry Smith   return 0;
60817ab2063SBarry Smith }
60917ab2063SBarry Smith 
610416022c9SBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,
61117ab2063SBarry Smith                              int *nzalloc,int *mem)
61217ab2063SBarry Smith {
613416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
614416022c9SBarry Smith   *nz      = a->nz;
615416022c9SBarry Smith   *nzalloc = a->maxnz;
616416022c9SBarry Smith   *mem     = (int)A->mem;
61717ab2063SBarry Smith   return 0;
61817ab2063SBarry Smith }
61917ab2063SBarry Smith 
62017ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
62117ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
62217ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
62317ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
62417ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
62517ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
62617ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
62717ab2063SBarry Smith 
62817ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
62917ab2063SBarry Smith {
630416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
631416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
63217ab2063SBarry Smith 
63317ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
63417ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
63517ab2063SBarry Smith   if (diag) {
63617ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
637416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
638416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
639416022c9SBarry Smith         a->ilen[rows[i]] = 1;
640416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
641416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
64217ab2063SBarry Smith       }
64317ab2063SBarry Smith       else {
64417ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
64517ab2063SBarry Smith         CHKERRQ(ierr);
64617ab2063SBarry Smith       }
64717ab2063SBarry Smith     }
64817ab2063SBarry Smith   }
64917ab2063SBarry Smith   else {
65017ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
651416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
652416022c9SBarry Smith       a->ilen[rows[i]] = 0;
65317ab2063SBarry Smith     }
65417ab2063SBarry Smith   }
65517ab2063SBarry Smith   ISRestoreIndices(is,&rows);
65617ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
65717ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
65817ab2063SBarry Smith   return 0;
65917ab2063SBarry Smith }
66017ab2063SBarry Smith 
661416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
66217ab2063SBarry Smith {
663416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
664416022c9SBarry Smith   *m = a->m; *n = a->n;
66517ab2063SBarry Smith   return 0;
66617ab2063SBarry Smith }
66717ab2063SBarry Smith 
668416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
66917ab2063SBarry Smith {
670416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
671416022c9SBarry Smith   *m = 0; *n = a->m;
67217ab2063SBarry Smith   return 0;
67317ab2063SBarry Smith }
674416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
67517ab2063SBarry Smith {
676416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
677416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
67817ab2063SBarry Smith 
679416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
68017ab2063SBarry Smith 
681416022c9SBarry Smith   if (!a->assembled) {
682416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
683416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
68417ab2063SBarry Smith   }
685416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
686416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
68717ab2063SBarry Smith   if (idx) {
68817ab2063SBarry Smith     if (*nz) {
689416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
69017ab2063SBarry Smith       *idx = (int *) PETSCMALLOC( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
69117ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
69217ab2063SBarry Smith     }
69317ab2063SBarry Smith     else *idx = 0;
69417ab2063SBarry Smith   }
69517ab2063SBarry Smith   return 0;
69617ab2063SBarry Smith }
69717ab2063SBarry Smith 
698416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
69917ab2063SBarry Smith {
70017ab2063SBarry Smith   if (idx) {if (*idx) PETSCFREE(*idx);}
70117ab2063SBarry Smith   return 0;
70217ab2063SBarry Smith }
70317ab2063SBarry Smith 
704416022c9SBarry Smith static int MatNorm_SeqAIJ(Mat A,MatNormType type,double *norm)
70517ab2063SBarry Smith {
706416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
707416022c9SBarry Smith   Scalar     *v = a->a;
70817ab2063SBarry Smith   double     sum = 0.0;
709416022c9SBarry Smith   int        i, j,shift = a->indexshift;
71017ab2063SBarry Smith 
711416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
71217ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
713416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
71417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
71517ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
71617ab2063SBarry Smith #else
71717ab2063SBarry Smith       sum += (*v)*(*v); v++;
71817ab2063SBarry Smith #endif
71917ab2063SBarry Smith     }
72017ab2063SBarry Smith     *norm = sqrt(sum);
72117ab2063SBarry Smith   }
72217ab2063SBarry Smith   else if (type == NORM_1) {
72317ab2063SBarry Smith     double *tmp;
724416022c9SBarry Smith     int    *jj = a->j;
725416022c9SBarry Smith     tmp = (double *) PETSCMALLOC( a->n*sizeof(double) ); CHKPTRQ(tmp);
726416022c9SBarry Smith     PetscZero(tmp,a->n*sizeof(double));
72717ab2063SBarry Smith     *norm = 0.0;
728416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
72917ab2063SBarry Smith #if defined(PETSC_COMPLEX)
73017ab2063SBarry Smith         tmp[*jj++ + shift] += abs(*v++);
73117ab2063SBarry Smith #else
73217ab2063SBarry Smith         tmp[*jj++ + shift] += fabs(*v++);
73317ab2063SBarry Smith #endif
73417ab2063SBarry Smith     }
735416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
73617ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
73717ab2063SBarry Smith     }
73817ab2063SBarry Smith     PETSCFREE(tmp);
73917ab2063SBarry Smith   }
74017ab2063SBarry Smith   else if (type == NORM_INFINITY) {
74117ab2063SBarry Smith     *norm = 0.0;
742416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
743416022c9SBarry Smith       v = a->a + a->i[j] + shift;
74417ab2063SBarry Smith       sum = 0.0;
745416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
74617ab2063SBarry Smith #if defined(PETSC_COMPLEX)
74717ab2063SBarry Smith         sum += abs(*v); v++;
74817ab2063SBarry Smith #else
74917ab2063SBarry Smith         sum += fabs(*v); v++;
75017ab2063SBarry Smith #endif
75117ab2063SBarry Smith       }
75217ab2063SBarry Smith       if (sum > *norm) *norm = sum;
75317ab2063SBarry Smith     }
75417ab2063SBarry Smith   }
75517ab2063SBarry Smith   else {
756*48d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
75717ab2063SBarry Smith   }
75817ab2063SBarry Smith   return 0;
75917ab2063SBarry Smith }
76017ab2063SBarry Smith 
761416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
76217ab2063SBarry Smith {
763416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
764416022c9SBarry Smith   Mat        C;
765416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
766416022c9SBarry Smith   Scalar     *array = a->a;
767416022c9SBarry Smith   int        shift = a->indexshift;
76817ab2063SBarry Smith 
769416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
770416022c9SBarry Smith   col = (int *) PETSCMALLOC((1+a->n)*sizeof(int)); CHKPTRQ(col);
771416022c9SBarry Smith   PetscZero(col,(1+a->n)*sizeof(int));
77217ab2063SBarry Smith   if (shift) {
77317ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
77417ab2063SBarry Smith   }
77517ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
776416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
77717ab2063SBarry Smith   PETSCFREE(col);
77817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
77917ab2063SBarry Smith     len = ai[i+1]-ai[i];
780416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
78117ab2063SBarry Smith     array += len; aj += len;
78217ab2063SBarry Smith   }
78317ab2063SBarry Smith   if (shift) {
78417ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
78517ab2063SBarry Smith   }
78617ab2063SBarry Smith 
787416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
788416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
78917ab2063SBarry Smith 
790416022c9SBarry Smith   if (B) {
791416022c9SBarry Smith     *B = C;
79217ab2063SBarry Smith   } else {
793416022c9SBarry Smith     /* This isn't really an in-place transpose */
794416022c9SBarry Smith     PETSCFREE(a->a);
795416022c9SBarry Smith     if (!a->singlemalloc) {PETSCFREE(a->i); PETSCFREE(a->j);}
796416022c9SBarry Smith     if (a->diag) PETSCFREE(a->diag);
797416022c9SBarry Smith     if (a->ilen) PETSCFREE(a->ilen);
798416022c9SBarry Smith     if (a->imax) PETSCFREE(a->imax);
799416022c9SBarry Smith     if (a->solve_work) PETSCFREE(a->solve_work);
800416022c9SBarry Smith     PETSCFREE(a);
801416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
802416022c9SBarry Smith     PETSCHEADERDESTROY(C);
80317ab2063SBarry Smith   }
80417ab2063SBarry Smith   return 0;
80517ab2063SBarry Smith }
80617ab2063SBarry Smith 
807416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
80817ab2063SBarry Smith {
809416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
81017ab2063SBarry Smith   Scalar     *l,*r,x,*v;
811416022c9SBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj;
812416022c9SBarry Smith   int        shift = a->indexshift;
81317ab2063SBarry Smith 
814*48d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
81517ab2063SBarry Smith   if (ll) {
81617ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
817416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
818416022c9SBarry Smith     v = a->a;
81917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
82017ab2063SBarry Smith       x = l[i];
821416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
82217ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
82317ab2063SBarry Smith     }
82417ab2063SBarry Smith   }
82517ab2063SBarry Smith   if (rr) {
82617ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
827416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
828416022c9SBarry Smith     v = a->a; jj = a->j;
82917ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
83017ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
83117ab2063SBarry Smith     }
83217ab2063SBarry Smith   }
83317ab2063SBarry Smith   return 0;
83417ab2063SBarry Smith }
83517ab2063SBarry Smith 
836416022c9SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,Mat *B)
83717ab2063SBarry Smith {
838416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
839416022c9SBarry Smith   int        nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n;
840416022c9SBarry Smith   int        *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift;
84117ab2063SBarry Smith   Scalar     *vwork;
842416022c9SBarry Smith   Mat        C;
84317ab2063SBarry Smith 
844416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
84517ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
84617ab2063SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
84717ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
84817ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
84917ab2063SBarry Smith 
85017ab2063SBarry Smith   smap  = (int *) PETSCMALLOC((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
85117ab2063SBarry Smith   cwork = (int *) PETSCMALLOC((1+ncols)*sizeof(int)); CHKPTRQ(cwork);
85217ab2063SBarry Smith   vwork = (Scalar *) PETSCMALLOC((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
853416022c9SBarry Smith   PetscZero(smap,oldcols*sizeof(int));
85417ab2063SBarry Smith   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
85517ab2063SBarry Smith 
85617ab2063SBarry Smith   /* Create and fill new matrix */
857416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,0,&C);CHKERRQ(ierr);
85817ab2063SBarry Smith   for (i=0; i<nrows; i++) {
85917ab2063SBarry Smith     nznew  = 0;
860416022c9SBarry Smith     kstart = a->i[irow[i]]+shift;
861416022c9SBarry Smith     kend   = kstart + a->ilen[irow[i]];
86217ab2063SBarry Smith     for ( k=kstart; k<kend; k++ ) {
863416022c9SBarry Smith       if (smap[a->j[k]+shift]) {
864416022c9SBarry Smith         cwork[nznew]   = smap[a->j[k]+shift] - 1;
865416022c9SBarry Smith         vwork[nznew++] = a->a[k];
86617ab2063SBarry Smith       }
86717ab2063SBarry Smith     }
868416022c9SBarry Smith     ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
86917ab2063SBarry Smith   }
870416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
871416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
87217ab2063SBarry Smith 
87317ab2063SBarry Smith   /* Free work space */
87417ab2063SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
87517ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
87617ab2063SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
877416022c9SBarry Smith   *B = C;
87817ab2063SBarry Smith   return 0;
87917ab2063SBarry Smith }
88017ab2063SBarry Smith 
88117ab2063SBarry Smith /* -------------------------------------------------------------------*/
88217ab2063SBarry Smith extern int MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,int,Mat *);
88317ab2063SBarry Smith extern int MatConvert_SeqAIJ(Mat,MatType,Mat *);
88417ab2063SBarry Smith static int MatCopyPrivate_SeqAIJ(Mat,Mat *);
88517ab2063SBarry Smith 
88617ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
88717ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
888416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
889416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
89017ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
89117ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
89217ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
89317ab2063SBarry Smith        MatRelax_SeqAIJ,
89417ab2063SBarry Smith        MatTranspose_SeqAIJ,
89517ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
89617ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
89717ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
89817ab2063SBarry Smith        MatCompress_SeqAIJ,
89917ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
90017ab2063SBarry Smith        MatGetReordering_SeqAIJ,
90117ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
90217ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
90317ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
90417ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
905416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
90617ab2063SBarry Smith        MatCopyPrivate_SeqAIJ};
90717ab2063SBarry Smith 
90817ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
90917ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
91017ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
91117ab2063SBarry Smith 
91217ab2063SBarry Smith /*@C
91317ab2063SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
91417ab2063SBarry Smith    (the default uniprocessor PETSc format).
91517ab2063SBarry Smith 
91617ab2063SBarry Smith    Input Parameters:
91717ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
91817ab2063SBarry Smith .  m - number of rows
91917ab2063SBarry Smith .  n - number of columns
92017ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
92117ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
92217ab2063SBarry Smith 
92317ab2063SBarry Smith    Output Parameter:
924416022c9SBarry Smith .  A - the matrix
92517ab2063SBarry Smith 
92617ab2063SBarry Smith    Notes:
92717ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
92817ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
92917ab2063SBarry Smith    storage.  That is, the stored row and column indices begin at
93017ab2063SBarry Smith    one, not zero.
93117ab2063SBarry Smith 
93217ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
93317ab2063SBarry Smith    Set both nz and nnz to zero for PETSc to control dynamic memory
93417ab2063SBarry Smith    allocation.
93517ab2063SBarry Smith 
93617ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse
93717ab2063SBarry Smith 
93817ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
93917ab2063SBarry Smith @*/
940416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
94117ab2063SBarry Smith {
942416022c9SBarry Smith   Mat        B;
943416022c9SBarry Smith   Mat_SeqAIJ *b;
94417ab2063SBarry Smith   int        i,len,ierr;
945416022c9SBarry Smith   *A      = 0;
946416022c9SBarry Smith   PETSCHEADERCREATE(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
947416022c9SBarry Smith   PLogObjectCreate(B);
948416022c9SBarry Smith   B->data               = (void *) (b = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(b);
949416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
950416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
951416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
952416022c9SBarry Smith   B->factor           = 0;
953416022c9SBarry Smith   B->lupivotthreshold = 1.0;
954416022c9SBarry Smith   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
955416022c9SBarry Smith   b->row              = 0;
956416022c9SBarry Smith   b->col              = 0;
957416022c9SBarry Smith   b->indexshift       = 0;
958416022c9SBarry Smith   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
95917ab2063SBarry Smith 
960416022c9SBarry Smith   b->m       = m;
961416022c9SBarry Smith   b->n       = n;
962416022c9SBarry Smith   b->imax    = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
96317ab2063SBarry Smith   if (!nnz) {
96417ab2063SBarry Smith     if (nz <= 0) nz = 1;
965416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
96617ab2063SBarry Smith     nz = nz*m;
96717ab2063SBarry Smith   }
96817ab2063SBarry Smith   else {
96917ab2063SBarry Smith     nz = 0;
970416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
97117ab2063SBarry Smith   }
97217ab2063SBarry Smith 
97317ab2063SBarry Smith   /* allocate the matrix space */
974416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
975416022c9SBarry Smith   b->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(b->a);
976416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
977416022c9SBarry Smith   PetscZero(b->j,nz*sizeof(int));
978416022c9SBarry Smith   b->i  = b->j + nz;
979416022c9SBarry Smith   b->singlemalloc = 1;
98017ab2063SBarry Smith 
981416022c9SBarry Smith   b->i[0] = -b->indexshift;
98217ab2063SBarry Smith   for (i=1; i<m+1; i++) {
983416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
98417ab2063SBarry Smith   }
98517ab2063SBarry Smith 
986416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
987416022c9SBarry Smith   b->ilen = (int *) PETSCMALLOC((m+1)*sizeof(int));
988416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
989416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
99017ab2063SBarry Smith 
991416022c9SBarry Smith   b->nz          = 0;
992416022c9SBarry Smith   b->maxnz       = nz;
993416022c9SBarry Smith   b->sorted      = 0;
994416022c9SBarry Smith   b->roworiented = 1;
995416022c9SBarry Smith   b->nonew       = 0;
996416022c9SBarry Smith   b->diag        = 0;
997416022c9SBarry Smith   b->assembled   = 0;
998416022c9SBarry Smith   b->solve_work  = 0;
99917ab2063SBarry Smith 
1000416022c9SBarry Smith   *A = B;
100117ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_superlu")) {
1002416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
100317ab2063SBarry Smith   }
100417ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_essl")) {
1005416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
100617ab2063SBarry Smith   }
100717ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_dxml")) {
1008416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1009416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
101017ab2063SBarry Smith   }
101117ab2063SBarry Smith 
101217ab2063SBarry Smith   return 0;
101317ab2063SBarry Smith }
101417ab2063SBarry Smith 
1015416022c9SBarry Smith static int MatCopyPrivate_SeqAIJ(Mat A,Mat *B)
101617ab2063SBarry Smith {
1017416022c9SBarry Smith   Mat        C;
1018416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1019416022c9SBarry Smith   int        i,len, m = a->m;
1020416022c9SBarry Smith   int        shift = a->indexshift;
1021416022c9SBarry Smith   *B      = 0;
102217ab2063SBarry Smith 
1023416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
1024416022c9SBarry Smith   PETSCHEADERCREATE(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1025416022c9SBarry Smith   PLogObjectCreate(C);
1026416022c9SBarry Smith   C->data       = (void *) (c = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(c);
1027416022c9SBarry Smith   PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps));
1028416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1029416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1030416022c9SBarry Smith   C->factor     = A->factor;
1031416022c9SBarry Smith   c->row        = 0;
1032416022c9SBarry Smith   c->col        = 0;
1033416022c9SBarry Smith   c->indexshift = shift;
103417ab2063SBarry Smith 
1035416022c9SBarry Smith   c->m          = a->m;
1036416022c9SBarry Smith   c->n          = a->n;
103717ab2063SBarry Smith 
1038416022c9SBarry Smith   c->imax       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1039416022c9SBarry Smith   c->ilen       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
104017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1041416022c9SBarry Smith     c->imax[i] = a->imax[i];
1042416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
104317ab2063SBarry Smith   }
104417ab2063SBarry Smith 
104517ab2063SBarry Smith   /* allocate the matrix space */
1046416022c9SBarry Smith   c->singlemalloc = 1;
1047416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1048416022c9SBarry Smith   c->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(c->a);
1049416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1050416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1051416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
105217ab2063SBarry Smith   if (m > 0) {
1053416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1054416022c9SBarry Smith     PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
105517ab2063SBarry Smith   }
105617ab2063SBarry Smith 
1057416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1058416022c9SBarry Smith   c->sorted      = a->sorted;
1059416022c9SBarry Smith   c->roworiented = a->roworiented;
1060416022c9SBarry Smith   c->nonew       = a->nonew;
106117ab2063SBarry Smith 
1062416022c9SBarry Smith   if (a->diag) {
1063416022c9SBarry Smith     c->diag = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1064416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
106517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1066416022c9SBarry Smith       c->diag[i] = a->diag[i];
106717ab2063SBarry Smith     }
106817ab2063SBarry Smith   }
1069416022c9SBarry Smith   else c->diag        = 0;
1070416022c9SBarry Smith   c->assembled        = 1;
1071416022c9SBarry Smith   c->nz               = a->nz;
1072416022c9SBarry Smith   c->maxnz            = a->maxnz;
1073416022c9SBarry Smith   c->solve_work       = 0;
1074416022c9SBarry Smith   *B = C;
107517ab2063SBarry Smith   return 0;
107617ab2063SBarry Smith }
107717ab2063SBarry Smith 
1078416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
107917ab2063SBarry Smith {
1080416022c9SBarry Smith   Mat_SeqAIJ   *a;
1081416022c9SBarry Smith   Mat          B;
1082416022c9SBarry Smith   int          i, nz, ierr, fd, header[4],numtid,*rowlengths = 0,M,N,shift;
108317ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
108417ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
108517ab2063SBarry Smith 
108617ab2063SBarry Smith   MPI_Comm_size(comm,&numtid);
108717ab2063SBarry Smith   if (numtid > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
108817ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1089416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
1090*48d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
109117ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
109217ab2063SBarry Smith 
109317ab2063SBarry Smith   /* read in row lengths */
109417ab2063SBarry Smith   rowlengths = (int*) PETSCMALLOC( M*sizeof(int) ); CHKPTRQ(rowlengths);
1095416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
109617ab2063SBarry Smith 
109717ab2063SBarry Smith   /* create our matrix */
1098416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1099416022c9SBarry Smith   B = *A;
1100416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1101416022c9SBarry Smith   shift = a->indexshift;
110217ab2063SBarry Smith 
110317ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1104416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
110517ab2063SBarry Smith   if (shift) {
110617ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1107416022c9SBarry Smith       a->j[i] += 1;
110817ab2063SBarry Smith     }
110917ab2063SBarry Smith   }
111017ab2063SBarry Smith 
111117ab2063SBarry Smith   /* read in nonzero values */
1112416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
111317ab2063SBarry Smith 
111417ab2063SBarry Smith   /* set matrix "i" values */
1115416022c9SBarry Smith   a->i[0] = -shift;
111617ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1117416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1118416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
111917ab2063SBarry Smith   }
112017ab2063SBarry Smith   PETSCFREE(rowlengths);
112117ab2063SBarry Smith 
1122416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1123416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
112417ab2063SBarry Smith   return 0;
112517ab2063SBarry Smith }
112617ab2063SBarry Smith 
112717ab2063SBarry Smith 
112817ab2063SBarry Smith 
1129