xref: /petsc/src/mat/impls/aij/seq/aij.c (revision e2f28af5970b8756dc9e817ba29989275be043b2)
117ab2063SBarry Smith #ifndef lint
2*e2f28af5SBarry Smith static char vcid[] = "$Id: aij.c,v 1.97 1995/10/11 22:09:34 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);
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 == COLUMNS_SORTED)            a->sorted      = 1;
330416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
331416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
332*e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
333*e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
334*e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
335*e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
336df719601SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
337df719601SLois Curfman McInnes   else if (op == COLUMN_ORIENTED)
338*e2f28af5SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED not supported");}
339df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
340df719601SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS not supported");}
341*e2f28af5SBarry Smith   else
342*e2f28af5SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:Option not supported");}
34317ab2063SBarry Smith   return 0;
34417ab2063SBarry Smith }
34517ab2063SBarry Smith 
346416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
34717ab2063SBarry Smith {
348416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
349416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
35017ab2063SBarry Smith   Scalar     *x, zero = 0.0;
35117ab2063SBarry Smith 
352416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
35317ab2063SBarry Smith   VecSet(&zero,v);
35417ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
355416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
356416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
357416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
358416022c9SBarry Smith       if (a->j[j]+shift == i) {
359416022c9SBarry Smith         x[i] = a->a[j];
36017ab2063SBarry Smith         break;
36117ab2063SBarry Smith       }
36217ab2063SBarry Smith     }
36317ab2063SBarry Smith   }
36417ab2063SBarry Smith   return 0;
36517ab2063SBarry Smith }
36617ab2063SBarry Smith 
36717ab2063SBarry Smith /* -------------------------------------------------------*/
36817ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
36917ab2063SBarry Smith /* -------------------------------------------------------*/
370416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
37117ab2063SBarry Smith {
372416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
37317ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
374416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
37517ab2063SBarry Smith 
376416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
37717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
378416022c9SBarry Smith   PetscZero(y,a->n*sizeof(Scalar));
37917ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
38017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
381416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
382416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
383416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
38417ab2063SBarry Smith     alpha = x[i];
38517ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
38617ab2063SBarry Smith   }
387416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
38817ab2063SBarry Smith   return 0;
38917ab2063SBarry Smith }
390416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
39117ab2063SBarry Smith {
392416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
39317ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
394416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
39517ab2063SBarry Smith 
396416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
39717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
39817ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
39917ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
40017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
401416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
402416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
403416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
40417ab2063SBarry Smith     alpha = x[i];
40517ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
40617ab2063SBarry Smith   }
40717ab2063SBarry Smith   return 0;
40817ab2063SBarry Smith }
40917ab2063SBarry Smith 
410416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
41117ab2063SBarry Smith {
412416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
41317ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
414416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
41517ab2063SBarry Smith 
416416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
41717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
41817ab2063SBarry Smith   x = x + shift; /* shift for Fortran start by 1 indexing */
419416022c9SBarry Smith   idx  = a->j;
420416022c9SBarry Smith   v    = a->a;
421416022c9SBarry Smith   ii   = a->i;
422416022c9SBarry Smith #if defined(PARCH_rs6000)
423416022c9SBarry Smith #pragma disjoint (*x,*v,*y)
424416022c9SBarry Smith #endif
42517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
426416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
42717ab2063SBarry Smith     sum  = 0.0;
42817ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
42917ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
430416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
43117ab2063SBarry Smith     y[i] = sum;
43217ab2063SBarry Smith   }
433416022c9SBarry Smith   PLogFlops(2*a->nz - m);
43417ab2063SBarry Smith   return 0;
43517ab2063SBarry Smith }
43617ab2063SBarry Smith 
437416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
43817ab2063SBarry Smith {
439416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
44017ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
441416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
44217ab2063SBarry Smith 
44348d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
44417ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
44517ab2063SBarry Smith   x = x + shift; /* shift for Fortran start by 1 indexing */
44617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
447416022c9SBarry Smith     idx  = a->j + a->i[i] + shift;
448416022c9SBarry Smith     v    = a->a + a->i[i] + shift;
449416022c9SBarry Smith     n    = a->i[i+1] - a->i[i];
45017ab2063SBarry Smith     sum  = y[i];
45117ab2063SBarry Smith     SPARSEDENSEDOT(sum,x,v,idx,n);
45217ab2063SBarry Smith     z[i] = sum;
45317ab2063SBarry Smith   }
454416022c9SBarry Smith   PLogFlops(2*a->nz);
45517ab2063SBarry Smith   return 0;
45617ab2063SBarry Smith }
45717ab2063SBarry Smith 
45817ab2063SBarry Smith /*
45917ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
46017ab2063SBarry Smith */
46117ab2063SBarry Smith 
462416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
46317ab2063SBarry Smith {
464416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
465416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
46617ab2063SBarry Smith 
467416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
46817ab2063SBarry Smith   diag = (int *) PETSCMALLOC( (m+1)*sizeof(int)); CHKPTRQ(diag);
469416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
470416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
471416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
472416022c9SBarry Smith       if (a->j[j]+shift == i) {
47317ab2063SBarry Smith         diag[i] = j - shift;
47417ab2063SBarry Smith         break;
47517ab2063SBarry Smith       }
47617ab2063SBarry Smith     }
47717ab2063SBarry Smith   }
478416022c9SBarry Smith   a->diag = diag;
47917ab2063SBarry Smith   return 0;
48017ab2063SBarry Smith }
48117ab2063SBarry Smith 
482416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
48317ab2063SBarry Smith                            double fshift,int its,Vec xx)
48417ab2063SBarry Smith {
485416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
486416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
487416022c9SBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i;
488416022c9SBarry Smith   int        shift = a->indexshift;
48917ab2063SBarry Smith 
49017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
491416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
492416022c9SBarry Smith   diag = a->diag;
493416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
49417ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
49517ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
49617ab2063SBarry Smith     bs = b + shift;
49717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
498416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
499416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
500416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
501416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
50217ab2063SBarry Smith         sum  = b[i]*d/omega;
50317ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
50417ab2063SBarry Smith         x[i] = sum;
50517ab2063SBarry Smith     }
50617ab2063SBarry Smith     return 0;
50717ab2063SBarry Smith   }
50817ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
50917ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
51017ab2063SBarry Smith   }
511416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
51217ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
51317ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
51417ab2063SBarry Smith 
51517ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
51617ab2063SBarry Smith 
51717ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
51817ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
51917ab2063SBarry Smith     is the relaxation factor.
52017ab2063SBarry Smith     */
52117ab2063SBarry Smith     t = (Scalar *) PETSCMALLOC( m*sizeof(Scalar) ); CHKPTRQ(t);
52217ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
52317ab2063SBarry Smith 
52417ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
52517ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
526416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
527416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
528416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
529416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
53017ab2063SBarry Smith       sum  = b[i];
53117ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
53217ab2063SBarry Smith       x[i] = omega*(sum/d);
53317ab2063SBarry Smith     }
53417ab2063SBarry Smith 
53517ab2063SBarry Smith     /*  t = b - (2*E - D)x */
536416022c9SBarry Smith     v = a->a;
53717ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
53817ab2063SBarry Smith 
53917ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
540416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
541416022c9SBarry Smith     diag = a->diag;
54217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
543416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
544416022c9SBarry Smith       n    = diag[i] - a->i[i];
545416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
546416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
54717ab2063SBarry Smith       sum  = t[i];
54817ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
54917ab2063SBarry Smith       t[i] = omega*(sum/d);
55017ab2063SBarry Smith     }
55117ab2063SBarry Smith 
55217ab2063SBarry Smith     /*  x = x + t */
55317ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
55417ab2063SBarry Smith     PETSCFREE(t);
55517ab2063SBarry Smith     return 0;
55617ab2063SBarry Smith   }
55717ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
55817ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
55917ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
560416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
561416022c9SBarry Smith         n    = diag[i] - a->i[i];
562416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
563416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
56417ab2063SBarry Smith         sum  = b[i];
56517ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
56617ab2063SBarry Smith         x[i] = omega*(sum/d);
56717ab2063SBarry Smith       }
56817ab2063SBarry Smith       xb = x;
56917ab2063SBarry Smith     }
57017ab2063SBarry Smith     else xb = b;
57117ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
57217ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
57317ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
574416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
57517ab2063SBarry Smith       }
57617ab2063SBarry Smith     }
57717ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
57817ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
579416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
580416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
581416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
582416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
58317ab2063SBarry Smith         sum  = xb[i];
58417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
58517ab2063SBarry Smith         x[i] = omega*(sum/d);
58617ab2063SBarry Smith       }
58717ab2063SBarry Smith     }
58817ab2063SBarry Smith     its--;
58917ab2063SBarry Smith   }
59017ab2063SBarry Smith   while (its--) {
59117ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
59217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
593416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
594416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
595416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
596416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
59717ab2063SBarry Smith         sum  = b[i];
59817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
59917ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
60017ab2063SBarry Smith       }
60117ab2063SBarry Smith     }
60217ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
60317ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
604416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
605416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
606416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
607416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
60817ab2063SBarry Smith         sum  = b[i];
60917ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
61017ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
61117ab2063SBarry Smith       }
61217ab2063SBarry Smith     }
61317ab2063SBarry Smith   }
61417ab2063SBarry Smith   return 0;
61517ab2063SBarry Smith }
61617ab2063SBarry Smith 
617416022c9SBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,
61817ab2063SBarry Smith                              int *nzalloc,int *mem)
61917ab2063SBarry Smith {
620416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
621416022c9SBarry Smith   *nz      = a->nz;
622416022c9SBarry Smith   *nzalloc = a->maxnz;
623416022c9SBarry Smith   *mem     = (int)A->mem;
62417ab2063SBarry Smith   return 0;
62517ab2063SBarry Smith }
62617ab2063SBarry Smith 
62717ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
62817ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
62917ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
63017ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
63117ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
63217ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
63317ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
63417ab2063SBarry Smith 
63517ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
63617ab2063SBarry Smith {
637416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
638416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
63917ab2063SBarry Smith 
64017ab2063SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
64117ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
64217ab2063SBarry Smith   if (diag) {
64317ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
644416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
645416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
646416022c9SBarry Smith         a->ilen[rows[i]] = 1;
647416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
648416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
64917ab2063SBarry Smith       }
65017ab2063SBarry Smith       else {
65117ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
65217ab2063SBarry Smith         CHKERRQ(ierr);
65317ab2063SBarry Smith       }
65417ab2063SBarry Smith     }
65517ab2063SBarry Smith   }
65617ab2063SBarry Smith   else {
65717ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
658416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
659416022c9SBarry Smith       a->ilen[rows[i]] = 0;
66017ab2063SBarry Smith     }
66117ab2063SBarry Smith   }
66217ab2063SBarry Smith   ISRestoreIndices(is,&rows);
66317ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
66417ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
66517ab2063SBarry Smith   return 0;
66617ab2063SBarry Smith }
66717ab2063SBarry Smith 
668416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
66917ab2063SBarry Smith {
670416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
671416022c9SBarry Smith   *m = a->m; *n = a->n;
67217ab2063SBarry Smith   return 0;
67317ab2063SBarry Smith }
67417ab2063SBarry Smith 
675416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
67617ab2063SBarry Smith {
677416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
678416022c9SBarry Smith   *m = 0; *n = a->m;
67917ab2063SBarry Smith   return 0;
68017ab2063SBarry Smith }
681416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
68217ab2063SBarry Smith {
683416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
684416022c9SBarry Smith   int        *itmp,i,ierr,shift = a->indexshift;
68517ab2063SBarry Smith 
686416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
68717ab2063SBarry Smith 
688416022c9SBarry Smith   if (!a->assembled) {
689416022c9SBarry Smith     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
690416022c9SBarry Smith     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
69117ab2063SBarry Smith   }
692416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
693416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
69417ab2063SBarry Smith   if (idx) {
69517ab2063SBarry Smith     if (*nz) {
696416022c9SBarry Smith       itmp = a->j + a->i[row] + shift;
69717ab2063SBarry Smith       *idx = (int *) PETSCMALLOC( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
69817ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
69917ab2063SBarry Smith     }
70017ab2063SBarry Smith     else *idx = 0;
70117ab2063SBarry Smith   }
70217ab2063SBarry Smith   return 0;
70317ab2063SBarry Smith }
70417ab2063SBarry Smith 
705416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
70617ab2063SBarry Smith {
70717ab2063SBarry Smith   if (idx) {if (*idx) PETSCFREE(*idx);}
70817ab2063SBarry Smith   return 0;
70917ab2063SBarry Smith }
71017ab2063SBarry Smith 
711416022c9SBarry Smith static int MatNorm_SeqAIJ(Mat A,MatNormType type,double *norm)
71217ab2063SBarry Smith {
713416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
714416022c9SBarry Smith   Scalar     *v = a->a;
71517ab2063SBarry Smith   double     sum = 0.0;
716416022c9SBarry Smith   int        i, j,shift = a->indexshift;
71717ab2063SBarry Smith 
718416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
71917ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
720416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
72117ab2063SBarry Smith #if defined(PETSC_COMPLEX)
72217ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
72317ab2063SBarry Smith #else
72417ab2063SBarry Smith       sum += (*v)*(*v); v++;
72517ab2063SBarry Smith #endif
72617ab2063SBarry Smith     }
72717ab2063SBarry Smith     *norm = sqrt(sum);
72817ab2063SBarry Smith   }
72917ab2063SBarry Smith   else if (type == NORM_1) {
73017ab2063SBarry Smith     double *tmp;
731416022c9SBarry Smith     int    *jj = a->j;
732416022c9SBarry Smith     tmp = (double *) PETSCMALLOC( a->n*sizeof(double) ); CHKPTRQ(tmp);
733416022c9SBarry Smith     PetscZero(tmp,a->n*sizeof(double));
73417ab2063SBarry Smith     *norm = 0.0;
735416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
73617ab2063SBarry Smith #if defined(PETSC_COMPLEX)
73717ab2063SBarry Smith         tmp[*jj++ + shift] += abs(*v++);
73817ab2063SBarry Smith #else
73917ab2063SBarry Smith         tmp[*jj++ + shift] += fabs(*v++);
74017ab2063SBarry Smith #endif
74117ab2063SBarry Smith     }
742416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
74317ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
74417ab2063SBarry Smith     }
74517ab2063SBarry Smith     PETSCFREE(tmp);
74617ab2063SBarry Smith   }
74717ab2063SBarry Smith   else if (type == NORM_INFINITY) {
74817ab2063SBarry Smith     *norm = 0.0;
749416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
750416022c9SBarry Smith       v = a->a + a->i[j] + shift;
75117ab2063SBarry Smith       sum = 0.0;
752416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
75317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
75417ab2063SBarry Smith         sum += abs(*v); v++;
75517ab2063SBarry Smith #else
75617ab2063SBarry Smith         sum += fabs(*v); v++;
75717ab2063SBarry Smith #endif
75817ab2063SBarry Smith       }
75917ab2063SBarry Smith       if (sum > *norm) *norm = sum;
76017ab2063SBarry Smith     }
76117ab2063SBarry Smith   }
76217ab2063SBarry Smith   else {
76348d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
76417ab2063SBarry Smith   }
76517ab2063SBarry Smith   return 0;
76617ab2063SBarry Smith }
76717ab2063SBarry Smith 
768416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
76917ab2063SBarry Smith {
770416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
771416022c9SBarry Smith   Mat        C;
772416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
773416022c9SBarry Smith   Scalar     *array = a->a;
774416022c9SBarry Smith   int        shift = a->indexshift;
77517ab2063SBarry Smith 
776416022c9SBarry Smith   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
777416022c9SBarry Smith   col = (int *) PETSCMALLOC((1+a->n)*sizeof(int)); CHKPTRQ(col);
778416022c9SBarry Smith   PetscZero(col,(1+a->n)*sizeof(int));
77917ab2063SBarry Smith   if (shift) {
78017ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
78117ab2063SBarry Smith   }
78217ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
783416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
78417ab2063SBarry Smith   PETSCFREE(col);
78517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
78617ab2063SBarry Smith     len = ai[i+1]-ai[i];
787416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
78817ab2063SBarry Smith     array += len; aj += len;
78917ab2063SBarry Smith   }
79017ab2063SBarry Smith   if (shift) {
79117ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
79217ab2063SBarry Smith   }
79317ab2063SBarry Smith 
794416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
795416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
79617ab2063SBarry Smith 
797416022c9SBarry Smith   if (B) {
798416022c9SBarry Smith     *B = C;
79917ab2063SBarry Smith   } else {
800416022c9SBarry Smith     /* This isn't really an in-place transpose */
801416022c9SBarry Smith     PETSCFREE(a->a);
802416022c9SBarry Smith     if (!a->singlemalloc) {PETSCFREE(a->i); PETSCFREE(a->j);}
803416022c9SBarry Smith     if (a->diag) PETSCFREE(a->diag);
804416022c9SBarry Smith     if (a->ilen) PETSCFREE(a->ilen);
805416022c9SBarry Smith     if (a->imax) PETSCFREE(a->imax);
806416022c9SBarry Smith     if (a->solve_work) PETSCFREE(a->solve_work);
807416022c9SBarry Smith     PETSCFREE(a);
808416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
809416022c9SBarry Smith     PETSCHEADERDESTROY(C);
81017ab2063SBarry Smith   }
81117ab2063SBarry Smith   return 0;
81217ab2063SBarry Smith }
81317ab2063SBarry Smith 
814416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
81517ab2063SBarry Smith {
816416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
81717ab2063SBarry Smith   Scalar     *l,*r,x,*v;
818416022c9SBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj;
819416022c9SBarry Smith   int        shift = a->indexshift;
82017ab2063SBarry Smith 
82148d91487SBarry Smith   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
82217ab2063SBarry Smith   if (ll) {
82317ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
824416022c9SBarry Smith     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
825416022c9SBarry Smith     v = a->a;
82617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
82717ab2063SBarry Smith       x = l[i];
828416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
82917ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
83017ab2063SBarry Smith     }
83117ab2063SBarry Smith   }
83217ab2063SBarry Smith   if (rr) {
83317ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
834416022c9SBarry Smith     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
835416022c9SBarry Smith     v = a->a; jj = a->j;
83617ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
83717ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
83817ab2063SBarry Smith     }
83917ab2063SBarry Smith   }
84017ab2063SBarry Smith   return 0;
84117ab2063SBarry Smith }
84217ab2063SBarry Smith 
843416022c9SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,Mat *B)
84417ab2063SBarry Smith {
845db02288aSLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c;
84602834360SBarry Smith   int        nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
84702834360SBarry Smith   int        *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
848db02288aSLois Curfman McInnes   int        first,step,*starts,*j_new,*i_new;
849db02288aSLois Curfman McInnes   Scalar     *vwork,*a_new;
850416022c9SBarry Smith   Mat        C;
85117ab2063SBarry Smith 
852416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
85317ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
85417ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
85517ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
85617ab2063SBarry Smith 
85702834360SBarry Smith   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
85802834360SBarry Smith     /* special case of contiguous rows */
859db02288aSLois Curfman McInnes     lens   = (int *) PETSCMALLOC((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
86002834360SBarry Smith     starts = lens + ncols;
86102834360SBarry Smith     /* loop over new rows determining lens and starting points */
86202834360SBarry Smith     for (i=0; i<nrows; i++) {
86302834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
86402834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
86502834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
86602834360SBarry Smith         if (a->j[k] >= first) {
86702834360SBarry Smith           starts[i] = k;
86802834360SBarry Smith           break;
86902834360SBarry Smith 	}
87002834360SBarry Smith       }
87102834360SBarry Smith       lens[i] = 0;
87202834360SBarry Smith       while (k < kend) {
87302834360SBarry Smith         if (a->j[k++] >= first+ncols) break;
87402834360SBarry Smith         lens[i]++;
87502834360SBarry Smith       }
87602834360SBarry Smith     }
87702834360SBarry Smith     /* create submatrix */
87802834360SBarry Smith     ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
879db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
880db02288aSLois Curfman McInnes 
88102834360SBarry Smith     /* loop over rows inserting into submatrix */
882db02288aSLois Curfman McInnes     a_new    = c->a;
883db02288aSLois Curfman McInnes     j_new    = c->j;
884db02288aSLois Curfman McInnes     i_new    = c->i;
885db02288aSLois Curfman McInnes     i_new[0] = -shift;
88602834360SBarry Smith     for (i=0; i<nrows; i++) {
88702834360SBarry Smith       for ( k=0; k<lens[i]; k++ ) {
888db02288aSLois Curfman McInnes         *j_new++ = a->j[k+starts[i]] - first;
88902834360SBarry Smith       }
890db02288aSLois Curfman McInnes       PetscMemcpy(a_new,a->a + starts[i],lens[i]*sizeof(Scalar));
891db02288aSLois Curfman McInnes       a_new      += lens[i];
892db02288aSLois Curfman McInnes       i_new[i+1]  = i_new[i] + lens[i];
8931987afe7SBarry Smith       c->ilen[i]  = lens[i];
89402834360SBarry Smith     }
895db02288aSLois Curfman McInnes     PETSCFREE(lens);
89602834360SBarry Smith   }
89702834360SBarry Smith   else {
89802834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
89917ab2063SBarry Smith     smap  = (int *) PETSCMALLOC((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
90002834360SBarry Smith     ssmap = smap + shift;
90102834360SBarry Smith     cwork = (int *) PETSCMALLOC((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork);
90202834360SBarry Smith     lens  = cwork + ncols;
90317ab2063SBarry Smith     vwork = (Scalar *) PETSCMALLOC((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
904416022c9SBarry Smith     PetscZero(smap,oldcols*sizeof(int));
90517ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
90602834360SBarry Smith     /* determine lens of each row */
90702834360SBarry Smith     for (i=0; i<nrows; i++) {
90802834360SBarry Smith       kstart  = a->i[irow[i]]+shift;
90902834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
91002834360SBarry Smith       lens[i] = 0;
91102834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
91202834360SBarry Smith         if (ssmap[a->j[k]]) {
91302834360SBarry Smith           lens[i]++;
91402834360SBarry Smith         }
91502834360SBarry Smith       }
91602834360SBarry Smith     }
91717ab2063SBarry Smith     /* Create and fill new matrix */
91802834360SBarry Smith     ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
91917ab2063SBarry Smith     for (i=0; i<nrows; i++) {
92017ab2063SBarry Smith       nznew  = 0;
921416022c9SBarry Smith       kstart = a->i[irow[i]]+shift;
922416022c9SBarry Smith       kend   = kstart + a->ilen[irow[i]];
92317ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
92402834360SBarry Smith         if (ssmap[a->j[k]]) {
92502834360SBarry Smith           cwork[nznew]   = ssmap[a->j[k]] - 1;
926416022c9SBarry Smith           vwork[nznew++] = a->a[k];
92717ab2063SBarry Smith         }
92817ab2063SBarry Smith       }
929416022c9SBarry Smith       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
93017ab2063SBarry Smith     }
93102834360SBarry Smith     /* Free work space */
93202834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
93302834360SBarry Smith     PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
93402834360SBarry Smith   }
935416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
936416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
93717ab2063SBarry Smith 
93817ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
939416022c9SBarry Smith   *B = C;
94017ab2063SBarry Smith   return 0;
94117ab2063SBarry Smith }
94217ab2063SBarry Smith 
94317ab2063SBarry Smith /* -------------------------------------------------------------------*/
94417ab2063SBarry Smith extern int MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,int,Mat *);
94517ab2063SBarry Smith extern int MatConvert_SeqAIJ(Mat,MatType,Mat *);
94617ab2063SBarry Smith static int MatCopyPrivate_SeqAIJ(Mat,Mat *);
94717ab2063SBarry Smith 
94817ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
94917ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
950416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
951416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
95217ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
95317ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
95417ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
95517ab2063SBarry Smith        MatRelax_SeqAIJ,
95617ab2063SBarry Smith        MatTranspose_SeqAIJ,
95717ab2063SBarry Smith        MatGetInfo_SeqAIJ,0,
95817ab2063SBarry Smith        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
95917ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
96017ab2063SBarry Smith        MatCompress_SeqAIJ,
96117ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
96217ab2063SBarry Smith        MatGetReordering_SeqAIJ,
96317ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
96417ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
96517ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
96617ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
967416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
96817ab2063SBarry Smith        MatCopyPrivate_SeqAIJ};
96917ab2063SBarry Smith 
97017ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
97117ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
97217ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
97317ab2063SBarry Smith 
97417ab2063SBarry Smith /*@C
97517ab2063SBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
97617ab2063SBarry Smith    (the default uniprocessor PETSc format).
97717ab2063SBarry Smith 
97817ab2063SBarry Smith    Input Parameters:
97917ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
98017ab2063SBarry Smith .  m - number of rows
98117ab2063SBarry Smith .  n - number of columns
98217ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
98317ab2063SBarry Smith .  nzz - number of nonzeros per row or null (possibly different for each row)
98417ab2063SBarry Smith 
98517ab2063SBarry Smith    Output Parameter:
986416022c9SBarry Smith .  A - the matrix
98717ab2063SBarry Smith 
98817ab2063SBarry Smith    Notes:
98917ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
99017ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
99117ab2063SBarry Smith    storage.  That is, the stored row and column indices begin at
99217ab2063SBarry Smith    one, not zero.
99317ab2063SBarry Smith 
99417ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
99517ab2063SBarry Smith    Set both nz and nnz to zero for PETSc to control dynamic memory
99617ab2063SBarry Smith    allocation.
99717ab2063SBarry Smith 
99817ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse
99917ab2063SBarry Smith 
100017ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
100117ab2063SBarry Smith @*/
1002416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
100317ab2063SBarry Smith {
1004416022c9SBarry Smith   Mat        B;
1005416022c9SBarry Smith   Mat_SeqAIJ *b;
100617ab2063SBarry Smith   int        i,len,ierr;
1007416022c9SBarry Smith   *A      = 0;
1008416022c9SBarry Smith   PETSCHEADERCREATE(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1009416022c9SBarry Smith   PLogObjectCreate(B);
1010416022c9SBarry Smith   B->data               = (void *) (b = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(b);
1011416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1012416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1013416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1014416022c9SBarry Smith   B->factor           = 0;
1015416022c9SBarry Smith   B->lupivotthreshold = 1.0;
1016416022c9SBarry Smith   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1017416022c9SBarry Smith   b->row              = 0;
1018416022c9SBarry Smith   b->col              = 0;
1019416022c9SBarry Smith   b->indexshift       = 0;
1020416022c9SBarry Smith   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
102117ab2063SBarry Smith 
1022416022c9SBarry Smith   b->m       = m;
1023416022c9SBarry Smith   b->n       = n;
1024416022c9SBarry Smith   b->imax    = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
102517ab2063SBarry Smith   if (!nnz) {
102617ab2063SBarry Smith     if (nz <= 0) nz = 1;
1027416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
102817ab2063SBarry Smith     nz = nz*m;
102917ab2063SBarry Smith   }
103017ab2063SBarry Smith   else {
103117ab2063SBarry Smith     nz = 0;
1032416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
103317ab2063SBarry Smith   }
103417ab2063SBarry Smith 
103517ab2063SBarry Smith   /* allocate the matrix space */
1036416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
1037416022c9SBarry Smith   b->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(b->a);
1038416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1039416022c9SBarry Smith   PetscZero(b->j,nz*sizeof(int));
1040416022c9SBarry Smith   b->i  = b->j + nz;
1041416022c9SBarry Smith   b->singlemalloc = 1;
104217ab2063SBarry Smith 
1043416022c9SBarry Smith   b->i[0] = -b->indexshift;
104417ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1045416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
104617ab2063SBarry Smith   }
104717ab2063SBarry Smith 
1048416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
1049416022c9SBarry Smith   b->ilen = (int *) PETSCMALLOC((m+1)*sizeof(int));
1050416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1051416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
105217ab2063SBarry Smith 
1053416022c9SBarry Smith   b->nz          = 0;
1054416022c9SBarry Smith   b->maxnz       = nz;
1055416022c9SBarry Smith   b->sorted      = 0;
1056416022c9SBarry Smith   b->roworiented = 1;
1057416022c9SBarry Smith   b->nonew       = 0;
1058416022c9SBarry Smith   b->diag        = 0;
1059416022c9SBarry Smith   b->assembled   = 0;
1060416022c9SBarry Smith   b->solve_work  = 0;
106117ab2063SBarry Smith 
1062416022c9SBarry Smith   *A = B;
106317ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_superlu")) {
1064416022c9SBarry Smith     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
106517ab2063SBarry Smith   }
106617ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_essl")) {
1067416022c9SBarry Smith     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
106817ab2063SBarry Smith   }
106917ab2063SBarry Smith   if (OptionsHasName(0,"-mat_aij_dxml")) {
1070416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1071416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
107217ab2063SBarry Smith   }
107317ab2063SBarry Smith 
107417ab2063SBarry Smith   return 0;
107517ab2063SBarry Smith }
107617ab2063SBarry Smith 
1077416022c9SBarry Smith static int MatCopyPrivate_SeqAIJ(Mat A,Mat *B)
107817ab2063SBarry Smith {
1079416022c9SBarry Smith   Mat        C;
1080416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1081416022c9SBarry Smith   int        i,len, m = a->m;
1082416022c9SBarry Smith   int        shift = a->indexshift;
1083416022c9SBarry Smith   *B      = 0;
108417ab2063SBarry Smith 
1085416022c9SBarry Smith   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
1086416022c9SBarry Smith   PETSCHEADERCREATE(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1087416022c9SBarry Smith   PLogObjectCreate(C);
1088416022c9SBarry Smith   C->data       = (void *) (c = PETSCNEW(Mat_SeqAIJ)); CHKPTRQ(c);
1089416022c9SBarry Smith   PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps));
1090416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1091416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1092416022c9SBarry Smith   C->factor     = A->factor;
1093416022c9SBarry Smith   c->row        = 0;
1094416022c9SBarry Smith   c->col        = 0;
1095416022c9SBarry Smith   c->indexshift = shift;
109617ab2063SBarry Smith 
1097416022c9SBarry Smith   c->m          = a->m;
1098416022c9SBarry Smith   c->n          = a->n;
109917ab2063SBarry Smith 
1100416022c9SBarry Smith   c->imax       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1101416022c9SBarry Smith   c->ilen       = (int *) PETSCMALLOC((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
110217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1103416022c9SBarry Smith     c->imax[i] = a->imax[i];
1104416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
110517ab2063SBarry Smith   }
110617ab2063SBarry Smith 
110717ab2063SBarry Smith   /* allocate the matrix space */
1108416022c9SBarry Smith   c->singlemalloc = 1;
1109416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1110416022c9SBarry Smith   c->a  = (Scalar *) PETSCMALLOC( len ); CHKPTRQ(c->a);
1111416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1112416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1113416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
111417ab2063SBarry Smith   if (m > 0) {
1115416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1116416022c9SBarry Smith     PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
111717ab2063SBarry Smith   }
111817ab2063SBarry Smith 
1119416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1120416022c9SBarry Smith   c->sorted      = a->sorted;
1121416022c9SBarry Smith   c->roworiented = a->roworiented;
1122416022c9SBarry Smith   c->nonew       = a->nonew;
112317ab2063SBarry Smith 
1124416022c9SBarry Smith   if (a->diag) {
1125416022c9SBarry Smith     c->diag = (int *) PETSCMALLOC( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1126416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
112717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1128416022c9SBarry Smith       c->diag[i] = a->diag[i];
112917ab2063SBarry Smith     }
113017ab2063SBarry Smith   }
1131416022c9SBarry Smith   else c->diag        = 0;
1132416022c9SBarry Smith   c->assembled        = 1;
1133416022c9SBarry Smith   c->nz               = a->nz;
1134416022c9SBarry Smith   c->maxnz            = a->maxnz;
1135416022c9SBarry Smith   c->solve_work       = 0;
1136416022c9SBarry Smith   *B = C;
113717ab2063SBarry Smith   return 0;
113817ab2063SBarry Smith }
113917ab2063SBarry Smith 
1140416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
114117ab2063SBarry Smith {
1142416022c9SBarry Smith   Mat_SeqAIJ   *a;
1143416022c9SBarry Smith   Mat          B;
1144416022c9SBarry Smith   int          i, nz, ierr, fd, header[4],numtid,*rowlengths = 0,M,N,shift;
114517ab2063SBarry Smith   PetscObject  vobj = (PetscObject) bview;
114617ab2063SBarry Smith   MPI_Comm     comm = vobj->comm;
114717ab2063SBarry Smith 
114817ab2063SBarry Smith   MPI_Comm_size(comm,&numtid);
114917ab2063SBarry Smith   if (numtid > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
115017ab2063SBarry Smith   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1151416022c9SBarry Smith   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
115248d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
115317ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
115417ab2063SBarry Smith 
115517ab2063SBarry Smith   /* read in row lengths */
115617ab2063SBarry Smith   rowlengths = (int*) PETSCMALLOC( M*sizeof(int) ); CHKPTRQ(rowlengths);
1157416022c9SBarry Smith   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
115817ab2063SBarry Smith 
115917ab2063SBarry Smith   /* create our matrix */
1160416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1161416022c9SBarry Smith   B = *A;
1162416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1163416022c9SBarry Smith   shift = a->indexshift;
116417ab2063SBarry Smith 
116517ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
1166416022c9SBarry Smith   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
116717ab2063SBarry Smith   if (shift) {
116817ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1169416022c9SBarry Smith       a->j[i] += 1;
117017ab2063SBarry Smith     }
117117ab2063SBarry Smith   }
117217ab2063SBarry Smith 
117317ab2063SBarry Smith   /* read in nonzero values */
1174416022c9SBarry Smith   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
117517ab2063SBarry Smith 
117617ab2063SBarry Smith   /* set matrix "i" values */
1177416022c9SBarry Smith   a->i[0] = -shift;
117817ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1179416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1180416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
118117ab2063SBarry Smith   }
118217ab2063SBarry Smith   PETSCFREE(rowlengths);
118317ab2063SBarry Smith 
1184416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1185416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
118617ab2063SBarry Smith   return 0;
118717ab2063SBarry Smith }
118817ab2063SBarry Smith 
118917ab2063SBarry Smith 
119017ab2063SBarry Smith 
1191