xref: /petsc/src/mat/impls/aij/seq/aij.c (revision a2ce50c755a9507c16aedae7d26067781b8c55ad)
16d84be18SBarry Smith 
217ab2063SBarry Smith #ifndef lint
3*a2ce50c7SBarry Smith static char vcid[] = "$Id: aij.c,v 1.173 1996/06/27 02:46:58 bsmith Exp bsmith $";
417ab2063SBarry Smith #endif
517ab2063SBarry Smith 
6d5d45c9bSBarry Smith /*
7d5d45c9bSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
8d5d45c9bSBarry Smith   matrix storage format.
9d5d45c9bSBarry Smith */
1017ab2063SBarry Smith #include "aij.h"
11f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
12f5eb4b81SSatish Balay #include "src/inline/spops.h"
13e4d965acSSatish Balay #include "petsc.h"
14f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1517ab2063SBarry Smith 
16*a2ce50c7SBarry Smith /*
17*a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
18*a2ce50c7SBarry Smith */
19*a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
20*a2ce50c7SBarry Smith {
21*a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
22*a2ce50c7SBarry Smith   int        ierr = 1;
23*a2ce50c7SBarry Smith 
24*a2ce50c7SBarry Smith   SETERRQ(ierr,"MatILUDTFactor_SeqAIJ:Not implemented");
25*a2ce50c7SBarry Smith }
26*a2ce50c7SBarry Smith 
27*a2ce50c7SBarry Smith /*
28*a2ce50c7SBarry Smith      Basic flow based reordering routine for AIJ storage format
29*a2ce50c7SBarry Smith */
30*a2ce50c7SBarry Smith int MatGetReordering_SeqAIJ_Flow(Mat A,IS *rperm,IS *cperm)
31*a2ce50c7SBarry Smith {
32*a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
33*a2ce50c7SBarry Smith   int        ierr = 1;
34*a2ce50c7SBarry Smith 
35*a2ce50c7SBarry Smith   SETERRQ(ierr,"MatGetReordering_SeqAIJ_Flow:Not implemented");
36*a2ce50c7SBarry Smith }
37*a2ce50c7SBarry Smith 
38bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3917ab2063SBarry Smith 
40*a2ce50c7SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatReordering type,IS *rperm, IS *cperm)
4117ab2063SBarry Smith {
42416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
43bcd2baecSBarry Smith   int        ierr, *ia, *ja,n,*idx,i,oshift,ishift;
4417ab2063SBarry Smith 
45a2744918SBarry Smith   /*
46*a2ce50c7SBarry Smith      This is tacky. The problem is we cannot use the registered ordering
47*a2ce50c7SBarry Smith     routines since they only work on nozero patterns. To have general
48*a2ce50c7SBarry Smith     registration means you register a particular ordering method for a
49*a2ce50c7SBarry Smith     particular data structure; hence something like a two dimensional
50*a2ce50c7SBarry Smith     look up table.
51*a2ce50c7SBarry Smith   */
52*a2ce50c7SBarry Smith   if (type == ORDER_FLOW) {
53*a2ce50c7SBarry Smith      return MatGetReordering_SeqAIJ_Flow(A,rperm,cperm);
54*a2ce50c7SBarry Smith   }
55*a2ce50c7SBarry Smith 
56*a2ce50c7SBarry Smith   /*
57a2744918SBarry Smith      this is tacky: In the future when we have written special factorization
58a2744918SBarry Smith      and solve routines for the identity permutation we should use a
59a2744918SBarry Smith      stride index set instead of the general one.
60a2744918SBarry Smith   */
61a2744918SBarry Smith   if (type  == ORDER_NATURAL) {
62a2744918SBarry Smith     n = a->n;
63a2744918SBarry Smith     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
64a2744918SBarry Smith     for ( i=0; i<n; i++ ) idx[i] = i;
65a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
66a2744918SBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
67a2744918SBarry Smith     PetscFree(idx);
68a2744918SBarry Smith     ISSetPermutation(*rperm);
69a2744918SBarry Smith     ISSetPermutation(*cperm);
70a2744918SBarry Smith     ISSetIdentity(*rperm);
71a2744918SBarry Smith     ISSetIdentity(*cperm);
72a2744918SBarry Smith     return 0;
73a2744918SBarry Smith   }
74a2744918SBarry Smith 
75bcd2baecSBarry Smith   MatReorderingRegisterAll();
76bcd2baecSBarry Smith   ishift = a->indexshift;
77bcd2baecSBarry Smith   oshift = -MatReorderingIndexShift[(int)type];
78bcd2baecSBarry Smith   if (MatReorderingRequiresSymmetric[(int)type]) {
79bcd2baecSBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja);
80bcd2baecSBarry Smith     CHKERRQ(ierr);
81416022c9SBarry Smith     ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
820452661fSBarry Smith     PetscFree(ia); PetscFree(ja);
83bcd2baecSBarry Smith   } else {
84bcd2baecSBarry Smith     if (ishift == oshift) {
85bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
86bcd2baecSBarry Smith     }
87bcd2baecSBarry Smith     else if (ishift == -1) {
88bcd2baecSBarry Smith       /* temporarily subtract 1 from i and j indices */
89bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
90bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
91bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
92bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
93bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
94bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
95bcd2baecSBarry Smith     } else {
96bcd2baecSBarry Smith       /* temporarily add 1 to i and j indices */
97bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
98bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
99bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
100bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr);
101bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
102bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
103bcd2baecSBarry Smith     }
104bcd2baecSBarry Smith   }
10517ab2063SBarry Smith   return 0;
10617ab2063SBarry Smith }
10717ab2063SBarry Smith 
108227d817aSBarry Smith #define CHUNKSIZE   15
10917ab2063SBarry Smith 
11017ab2063SBarry Smith /* This version has row oriented v  */
111416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
11217ab2063SBarry Smith {
113416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
114416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1154b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
116d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
117416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
11817ab2063SBarry Smith 
11917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
120416022c9SBarry Smith     row  = im[k];
12117ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
122416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
12317ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
12417ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
125416022c9SBarry Smith     low = 0;
12617ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
127416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
128416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
1294b0e389bSBarry Smith       col = in[l] - shift;
1304b0e389bSBarry Smith       if (roworiented) {
1314b0e389bSBarry Smith         value = *v++;
1324b0e389bSBarry Smith       }
1334b0e389bSBarry Smith       else {
1344b0e389bSBarry Smith         value = v[k + l*m];
1354b0e389bSBarry Smith       }
136416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
137416022c9SBarry Smith       while (high-low > 5) {
138416022c9SBarry Smith         t = (low+high)/2;
139416022c9SBarry Smith         if (rp[t] > col) high = t;
140416022c9SBarry Smith         else             low  = t;
14117ab2063SBarry Smith       }
142416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
14317ab2063SBarry Smith         if (rp[i] > col) break;
14417ab2063SBarry Smith         if (rp[i] == col) {
145416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
14617ab2063SBarry Smith           else                  ap[i] = value;
14717ab2063SBarry Smith           goto noinsert;
14817ab2063SBarry Smith         }
14917ab2063SBarry Smith       }
15017ab2063SBarry Smith       if (nonew) goto noinsert;
15117ab2063SBarry Smith       if (nrow >= rmax) {
15217ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
153416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
15417ab2063SBarry Smith         Scalar *new_a;
15517ab2063SBarry Smith 
15617ab2063SBarry Smith         /* malloc new storage space */
157416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1580452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
15917ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
16017ab2063SBarry Smith         new_i   = new_j + new_nz;
16117ab2063SBarry Smith 
16217ab2063SBarry Smith         /* copy over old data into new slots */
16317ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
164416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
165416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
166416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
167416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
16817ab2063SBarry Smith                                                            len*sizeof(int));
169416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
170416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
17117ab2063SBarry Smith                                                            len*sizeof(Scalar));
17217ab2063SBarry Smith         /* free up old matrix storage */
1730452661fSBarry Smith         PetscFree(a->a);
1740452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
175416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
176416022c9SBarry Smith         a->singlemalloc = 1;
17717ab2063SBarry Smith 
17817ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
179416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
180416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
181416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
182b810aeb4SBarry Smith         a->reallocs++;
18317ab2063SBarry Smith       }
184416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
185416022c9SBarry Smith       /* shift up all the later entries in this row */
186416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
18717ab2063SBarry Smith         rp[ii+1] = rp[ii];
18817ab2063SBarry Smith         ap[ii+1] = ap[ii];
18917ab2063SBarry Smith       }
19017ab2063SBarry Smith       rp[i] = col;
19117ab2063SBarry Smith       ap[i] = value;
19217ab2063SBarry Smith       noinsert:;
193416022c9SBarry Smith       low = i + 1;
19417ab2063SBarry Smith     }
19517ab2063SBarry Smith     ailen[row] = nrow;
19617ab2063SBarry Smith   }
19717ab2063SBarry Smith   return 0;
19817ab2063SBarry Smith }
19917ab2063SBarry Smith 
2007eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2017eb43aa7SLois Curfman McInnes {
2027eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
203b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2047eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2057eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2067eb43aa7SLois Curfman McInnes 
2077eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2087eb43aa7SLois Curfman McInnes     row  = im[k];
2097eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
2107eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
2117eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2127eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2137eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
2147eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
2157eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
2167eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2177eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2187eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2197eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2207eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2217eb43aa7SLois Curfman McInnes         else             low  = t;
2227eb43aa7SLois Curfman McInnes       }
2237eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2247eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2257eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
226b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2277eb43aa7SLois Curfman McInnes           goto finished;
2287eb43aa7SLois Curfman McInnes         }
2297eb43aa7SLois Curfman McInnes       }
230b49de8d1SLois Curfman McInnes       *v++ = zero;
2317eb43aa7SLois Curfman McInnes       finished:;
2327eb43aa7SLois Curfman McInnes     }
2337eb43aa7SLois Curfman McInnes   }
2347eb43aa7SLois Curfman McInnes   return 0;
2357eb43aa7SLois Curfman McInnes }
2367eb43aa7SLois Curfman McInnes 
23717ab2063SBarry Smith #include "draw.h"
23817ab2063SBarry Smith #include "pinclude/pviewer.h"
23977c4ece6SBarry Smith #include "sys.h"
24017ab2063SBarry Smith 
241416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
24217ab2063SBarry Smith {
243416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
244416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
24517ab2063SBarry Smith 
24690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2470452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
248416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
249416022c9SBarry Smith   col_lens[1] = a->m;
250416022c9SBarry Smith   col_lens[2] = a->n;
251416022c9SBarry Smith   col_lens[3] = a->nz;
252416022c9SBarry Smith 
253416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
254416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
255416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
25617ab2063SBarry Smith   }
25777c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
2580452661fSBarry Smith   PetscFree(col_lens);
259416022c9SBarry Smith 
260416022c9SBarry Smith   /* store column indices (zero start index) */
261416022c9SBarry Smith   if (a->indexshift) {
262416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
26317ab2063SBarry Smith   }
26477c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr);
265416022c9SBarry Smith   if (a->indexshift) {
266416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
26717ab2063SBarry Smith   }
268416022c9SBarry Smith 
269416022c9SBarry Smith   /* store nonzero values */
27077c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
27117ab2063SBarry Smith   return 0;
27217ab2063SBarry Smith }
273416022c9SBarry Smith 
274416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
275416022c9SBarry Smith {
276416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
277496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
27817ab2063SBarry Smith   FILE        *fd;
27917ab2063SBarry Smith   char        *outputname;
28017ab2063SBarry Smith 
28190ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
282416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
28390ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
28490ace30eSBarry Smith   if (format == ASCII_FORMAT_INFO) {
28595e01e2fSLois Curfman McInnes     return 0;
28695e01e2fSLois Curfman McInnes   }
28790ace30eSBarry Smith   else if (format == ASCII_FORMAT_INFO_DETAILED) {
288496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
289496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
290496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
29195e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
29295e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
29317ab2063SBarry Smith   }
29490ace30eSBarry Smith   else if (format == ASCII_FORMAT_MATLAB) {
29517ab2063SBarry Smith     int nz, nzalloc, mem;
296416022c9SBarry Smith     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
297416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
29817ab2063SBarry Smith     fprintf(fd,"%% Nonzeros = %d \n",nz);
29917ab2063SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
30017ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
30117ab2063SBarry Smith 
30217ab2063SBarry Smith     for (i=0; i<m; i++) {
303416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
30417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
3057a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]),
306416022c9SBarry Smith                    imag(a->a[j]));
30717ab2063SBarry Smith #else
3087a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
30917ab2063SBarry Smith #endif
31017ab2063SBarry Smith       }
31117ab2063SBarry Smith     }
31217ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
31317ab2063SBarry Smith   }
31444cd7ae7SLois Curfman McInnes   else if (format == ASCII_FORMAT_COMMON) {
31544cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
31644cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
31744cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
31844cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
31944cd7ae7SLois Curfman McInnes         if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0)
32044cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
32144cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
32244cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
32344cd7ae7SLois Curfman McInnes #else
32444cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
32544cd7ae7SLois Curfman McInnes #endif
32644cd7ae7SLois Curfman McInnes       }
32744cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
32844cd7ae7SLois Curfman McInnes     }
32944cd7ae7SLois Curfman McInnes   }
33017ab2063SBarry Smith   else {
33117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
33217ab2063SBarry Smith       fprintf(fd,"row %d:",i);
333416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
33417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
335416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
336416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
33717ab2063SBarry Smith         }
33817ab2063SBarry Smith         else {
339416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
34017ab2063SBarry Smith         }
34117ab2063SBarry Smith #else
342416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
34317ab2063SBarry Smith #endif
34417ab2063SBarry Smith       }
34517ab2063SBarry Smith       fprintf(fd,"\n");
34617ab2063SBarry Smith     }
34717ab2063SBarry Smith   }
34817ab2063SBarry Smith   fflush(fd);
349416022c9SBarry Smith   return 0;
350416022c9SBarry Smith }
351416022c9SBarry Smith 
352416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
353416022c9SBarry Smith {
354416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
355cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
356cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
357bcd2baecSBarry Smith   Draw        draw;
358cddf8d76SBarry Smith   DrawButton  button;
35919bcc07fSBarry Smith   PetscTruth  isnull;
360cddf8d76SBarry Smith 
361bcd2baecSBarry Smith   ViewerDrawGetDraw(viewer,&draw);
36219bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
36319bcc07fSBarry Smith 
364416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
365416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
366416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
367416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
368cddf8d76SBarry Smith   color = DRAW_BLUE;
369416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
370cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
371416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
372cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
373cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
374cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
375cddf8d76SBarry Smith #else
376cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
377cddf8d76SBarry Smith #endif
378cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
379cddf8d76SBarry Smith     }
380cddf8d76SBarry Smith   }
381cddf8d76SBarry Smith   color = DRAW_CYAN;
382cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
383cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
384cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
385cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
386cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
387cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
388cddf8d76SBarry Smith     }
389cddf8d76SBarry Smith   }
390cddf8d76SBarry Smith   color = DRAW_RED;
391cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
392cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
393cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
394cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
395cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
396cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
397cddf8d76SBarry Smith #else
398cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
399cddf8d76SBarry Smith #endif
400cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
401416022c9SBarry Smith     }
402416022c9SBarry Smith   }
403416022c9SBarry Smith   DrawFlush(draw);
404cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
405cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
406cddf8d76SBarry Smith 
407cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
408cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
409cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
410cddf8d76SBarry Smith     DrawClear(draw);
411cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
412cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
413cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
414cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
415cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
416cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
417cddf8d76SBarry Smith     w *= scale; h *= scale;
418cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
419cddf8d76SBarry Smith     color = DRAW_BLUE;
420cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
421cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
422cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
423cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
424cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
425cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
426cddf8d76SBarry Smith #else
427cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
428cddf8d76SBarry Smith #endif
429cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
430cddf8d76SBarry Smith       }
431cddf8d76SBarry Smith     }
432cddf8d76SBarry Smith     color = DRAW_CYAN;
433cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
434cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
435cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
436cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
437cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
438cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
439cddf8d76SBarry Smith       }
440cddf8d76SBarry Smith     }
441cddf8d76SBarry Smith     color = DRAW_RED;
442cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
443cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
444cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
445cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
446cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
447cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
448cddf8d76SBarry Smith #else
449cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
450cddf8d76SBarry Smith #endif
451cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
452cddf8d76SBarry Smith       }
453cddf8d76SBarry Smith     }
454cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
455cddf8d76SBarry Smith   }
456416022c9SBarry Smith   return 0;
457416022c9SBarry Smith }
458416022c9SBarry Smith 
459416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
460416022c9SBarry Smith {
461416022c9SBarry Smith   Mat         A = (Mat) obj;
462416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
463bcd2baecSBarry Smith   ViewerType  vtype;
464bcd2baecSBarry Smith   int         ierr;
465416022c9SBarry Smith 
466416022c9SBarry Smith   if (!viewer) {
467bcd2baecSBarry Smith     viewer = STDOUT_VIEWER_SELF;
468416022c9SBarry Smith   }
469bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
470bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
471416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
472416022c9SBarry Smith   }
473bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
474416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
475416022c9SBarry Smith   }
476bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
477416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
478416022c9SBarry Smith   }
479bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
480bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
48117ab2063SBarry Smith   }
48217ab2063SBarry Smith   return 0;
48317ab2063SBarry Smith }
48419bcc07fSBarry Smith 
485c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
486416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
48717ab2063SBarry Smith {
488416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
48941c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
490416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
491416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
49217ab2063SBarry Smith 
49317ab2063SBarry Smith   if (mode == FLUSH_ASSEMBLY) return 0;
49417ab2063SBarry Smith 
49517ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
496416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
49717ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
49817ab2063SBarry Smith     if (fshift) {
499416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
50017ab2063SBarry Smith       N = ailen[i];
50117ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
50217ab2063SBarry Smith         ip[j-fshift] = ip[j];
50317ab2063SBarry Smith         ap[j-fshift] = ap[j];
50417ab2063SBarry Smith       }
50517ab2063SBarry Smith     }
50617ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
50717ab2063SBarry Smith   }
50817ab2063SBarry Smith   if (m) {
50917ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
51017ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
51117ab2063SBarry Smith   }
51217ab2063SBarry Smith   /* reset ilen and imax for each row */
51317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
51417ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
51517ab2063SBarry Smith   }
516416022c9SBarry Smith   a->nz = ai[m] + shift;
51717ab2063SBarry Smith 
51817ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
519416022c9SBarry Smith   if (fshift && a->diag) {
5200452661fSBarry Smith     PetscFree(a->diag);
521416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
522416022c9SBarry Smith     a->diag = 0;
52317ab2063SBarry Smith   }
52494a424c1SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n",
525b810aeb4SBarry Smith            fshift,a->nz,m);
52694a424c1SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n",
527b810aeb4SBarry Smith            a->reallocs);
52876dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
52941c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
53017ab2063SBarry Smith   return 0;
53117ab2063SBarry Smith }
53217ab2063SBarry Smith 
533416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
53417ab2063SBarry Smith {
535416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
536cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
53717ab2063SBarry Smith   return 0;
53817ab2063SBarry Smith }
539416022c9SBarry Smith 
54017ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
54117ab2063SBarry Smith {
542416022c9SBarry Smith   Mat        A  = (Mat) obj;
543416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
544d5d45c9bSBarry Smith 
54517ab2063SBarry Smith #if defined(PETSC_LOG)
546416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
54717ab2063SBarry Smith #endif
5480452661fSBarry Smith   PetscFree(a->a);
5490452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
5500452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
5510452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
5520452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
5530452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
55476dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
5550452661fSBarry Smith   PetscFree(a);
556f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
557f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
55817ab2063SBarry Smith   return 0;
55917ab2063SBarry Smith }
56017ab2063SBarry Smith 
561416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
56217ab2063SBarry Smith {
56317ab2063SBarry Smith   return 0;
56417ab2063SBarry Smith }
56517ab2063SBarry Smith 
566416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
56717ab2063SBarry Smith {
568416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
569416022c9SBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
5704b0e389bSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
571416022c9SBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
572416022c9SBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
573416022c9SBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
574e2f28af5SBarry Smith   else if (op == ROWS_SORTED ||
575e2f28af5SBarry Smith            op == SYMMETRIC_MATRIX ||
576e2f28af5SBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
577e2f28af5SBarry Smith            op == YES_NEW_DIAGONALS)
57894a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
579df719601SLois Curfman McInnes   else if (op == NO_NEW_DIAGONALS)
5804d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
5816c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_1)            a->inode.limit  = 1;
5826c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_2)            a->inode.limit  = 2;
5836c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_3)            a->inode.limit  = 3;
5846c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_4)            a->inode.limit  = 4;
5856c7ebb05SLois Curfman McInnes   else if (op == INODE_LIMIT_5)            a->inode.limit  = 5;
586e2f28af5SBarry Smith   else
5874d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
58817ab2063SBarry Smith   return 0;
58917ab2063SBarry Smith }
59017ab2063SBarry Smith 
591416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
59217ab2063SBarry Smith {
593416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
594416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
59517ab2063SBarry Smith   Scalar     *x, zero = 0.0;
59617ab2063SBarry Smith 
59717ab2063SBarry Smith   VecSet(&zero,v);
59817ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
599416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
600416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
601416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
602416022c9SBarry Smith       if (a->j[j]+shift == i) {
603416022c9SBarry Smith         x[i] = a->a[j];
60417ab2063SBarry Smith         break;
60517ab2063SBarry Smith       }
60617ab2063SBarry Smith     }
60717ab2063SBarry Smith   }
60817ab2063SBarry Smith   return 0;
60917ab2063SBarry Smith }
61017ab2063SBarry Smith 
61117ab2063SBarry Smith /* -------------------------------------------------------*/
61217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
61317ab2063SBarry Smith /* -------------------------------------------------------*/
61444cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
61517ab2063SBarry Smith {
616416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
61717ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
618416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
61917ab2063SBarry Smith 
62017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
621cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
62217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
62317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
624416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
625416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
626416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
62717ab2063SBarry Smith     alpha = x[i];
62817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
62917ab2063SBarry Smith   }
630416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
63117ab2063SBarry Smith   return 0;
63217ab2063SBarry Smith }
633d5d45c9bSBarry Smith 
63444cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
63517ab2063SBarry Smith {
636416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
63717ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
638416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
63917ab2063SBarry Smith 
64017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
64117ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
64217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
64317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
644416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
645416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
646416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
64717ab2063SBarry Smith     alpha = x[i];
64817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
64917ab2063SBarry Smith   }
65017ab2063SBarry Smith   return 0;
65117ab2063SBarry Smith }
65217ab2063SBarry Smith 
65344cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
65417ab2063SBarry Smith {
655416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
65617ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
657416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
65817ab2063SBarry Smith 
65917ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
66017ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
661416022c9SBarry Smith   idx  = a->j;
662416022c9SBarry Smith   v    = a->a;
663416022c9SBarry Smith   ii   = a->i;
66417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
665416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
66617ab2063SBarry Smith     sum  = 0.0;
66717ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
66817ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
669416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
67017ab2063SBarry Smith     y[i] = sum;
67117ab2063SBarry Smith   }
672416022c9SBarry Smith   PLogFlops(2*a->nz - m);
67317ab2063SBarry Smith   return 0;
67417ab2063SBarry Smith }
67517ab2063SBarry Smith 
67644cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
67717ab2063SBarry Smith {
678416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
67917ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
680cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
68117ab2063SBarry Smith 
68217ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
68317ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
684cddf8d76SBarry Smith   idx  = a->j;
685cddf8d76SBarry Smith   v    = a->a;
686cddf8d76SBarry Smith   ii   = a->i;
68717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
688cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
68917ab2063SBarry Smith     sum  = y[i];
690cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
691cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
69217ab2063SBarry Smith     z[i] = sum;
69317ab2063SBarry Smith   }
694416022c9SBarry Smith   PLogFlops(2*a->nz);
69517ab2063SBarry Smith   return 0;
69617ab2063SBarry Smith }
69717ab2063SBarry Smith 
69817ab2063SBarry Smith /*
69917ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
70017ab2063SBarry Smith */
70117ab2063SBarry Smith 
702416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
70317ab2063SBarry Smith {
704416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
705416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
70617ab2063SBarry Smith 
7070452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
708416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
709416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
710416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
711416022c9SBarry Smith       if (a->j[j]+shift == i) {
71217ab2063SBarry Smith         diag[i] = j - shift;
71317ab2063SBarry Smith         break;
71417ab2063SBarry Smith       }
71517ab2063SBarry Smith     }
71617ab2063SBarry Smith   }
717416022c9SBarry Smith   a->diag = diag;
71817ab2063SBarry Smith   return 0;
71917ab2063SBarry Smith }
72017ab2063SBarry Smith 
72144cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
72217ab2063SBarry Smith                            double fshift,int its,Vec xx)
72317ab2063SBarry Smith {
724416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
725416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
726d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
72717ab2063SBarry Smith 
72817ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
729416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
730416022c9SBarry Smith   diag = a->diag;
731416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
73217ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
73317ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
73417ab2063SBarry Smith     bs = b + shift;
73517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
736416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
737416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
738416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
739416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
74017ab2063SBarry Smith         sum  = b[i]*d/omega;
74117ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
74217ab2063SBarry Smith         x[i] = sum;
74317ab2063SBarry Smith     }
74417ab2063SBarry Smith     return 0;
74517ab2063SBarry Smith   }
74617ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
74717ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
74817ab2063SBarry Smith   }
749416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
75017ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
75117ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
75217ab2063SBarry Smith 
75317ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
75417ab2063SBarry Smith 
75517ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
75617ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
75717ab2063SBarry Smith     is the relaxation factor.
75817ab2063SBarry Smith     */
7590452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
76017ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
76117ab2063SBarry Smith 
76217ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
76317ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
764416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
765416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
766416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
767416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
76817ab2063SBarry Smith       sum  = b[i];
76917ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
77017ab2063SBarry Smith       x[i] = omega*(sum/d);
77117ab2063SBarry Smith     }
77217ab2063SBarry Smith 
77317ab2063SBarry Smith     /*  t = b - (2*E - D)x */
774416022c9SBarry Smith     v = a->a;
77517ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
77617ab2063SBarry Smith 
77717ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
778416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
779416022c9SBarry Smith     diag = a->diag;
78017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
781416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
782416022c9SBarry Smith       n    = diag[i] - a->i[i];
783416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
784416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
78517ab2063SBarry Smith       sum  = t[i];
78617ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
78717ab2063SBarry Smith       t[i] = omega*(sum/d);
78817ab2063SBarry Smith     }
78917ab2063SBarry Smith 
79017ab2063SBarry Smith     /*  x = x + t */
79117ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7920452661fSBarry Smith     PetscFree(t);
79317ab2063SBarry Smith     return 0;
79417ab2063SBarry Smith   }
79517ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
79617ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
79717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
798416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
799416022c9SBarry Smith         n    = diag[i] - a->i[i];
800416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
801416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
80217ab2063SBarry Smith         sum  = b[i];
80317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
80417ab2063SBarry Smith         x[i] = omega*(sum/d);
80517ab2063SBarry Smith       }
80617ab2063SBarry Smith       xb = x;
80717ab2063SBarry Smith     }
80817ab2063SBarry Smith     else xb = b;
80917ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
81017ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
81117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
812416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
81317ab2063SBarry Smith       }
81417ab2063SBarry Smith     }
81517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
81617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
817416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
818416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
819416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
820416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
82117ab2063SBarry Smith         sum  = xb[i];
82217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
82317ab2063SBarry Smith         x[i] = omega*(sum/d);
82417ab2063SBarry Smith       }
82517ab2063SBarry Smith     }
82617ab2063SBarry Smith     its--;
82717ab2063SBarry Smith   }
82817ab2063SBarry Smith   while (its--) {
82917ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
83017ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
831416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
832416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
833416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
834416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
83517ab2063SBarry Smith         sum  = b[i];
83617ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
83717ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
83817ab2063SBarry Smith       }
83917ab2063SBarry Smith     }
84017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
84117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
842416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
843416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
844416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
845416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
84617ab2063SBarry Smith         sum  = b[i];
84717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
84817ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
84917ab2063SBarry Smith       }
85017ab2063SBarry Smith     }
85117ab2063SBarry Smith   }
85217ab2063SBarry Smith   return 0;
85317ab2063SBarry Smith }
85417ab2063SBarry Smith 
855d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
85617ab2063SBarry Smith {
857416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
858bcd2baecSBarry Smith   if (nz)      *nz      = a->nz;
859bcd2baecSBarry Smith   if (nzalloc) *nzalloc = a->maxnz;
860bcd2baecSBarry Smith   if (mem)     *mem     = (int)A->mem;
86117ab2063SBarry Smith   return 0;
86217ab2063SBarry Smith }
86317ab2063SBarry Smith 
86417ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
86517ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
86617ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
86717ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
86817ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
86917ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
87017ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
87117ab2063SBarry Smith 
87217ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
87317ab2063SBarry Smith {
874416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
875416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
87617ab2063SBarry Smith 
87777c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
87817ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
87917ab2063SBarry Smith   if (diag) {
88017ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
881416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
882416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
883416022c9SBarry Smith         a->ilen[rows[i]] = 1;
884416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
885416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
88617ab2063SBarry Smith       }
88717ab2063SBarry Smith       else {
88817ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
88917ab2063SBarry Smith         CHKERRQ(ierr);
89017ab2063SBarry Smith       }
89117ab2063SBarry Smith     }
89217ab2063SBarry Smith   }
89317ab2063SBarry Smith   else {
89417ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
895416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
896416022c9SBarry Smith       a->ilen[rows[i]] = 0;
89717ab2063SBarry Smith     }
89817ab2063SBarry Smith   }
89917ab2063SBarry Smith   ISRestoreIndices(is,&rows);
90017ab2063SBarry Smith   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
90117ab2063SBarry Smith   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
90217ab2063SBarry Smith   return 0;
90317ab2063SBarry Smith }
90417ab2063SBarry Smith 
905416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
90617ab2063SBarry Smith {
907416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
908416022c9SBarry Smith   *m = a->m; *n = a->n;
90917ab2063SBarry Smith   return 0;
91017ab2063SBarry Smith }
91117ab2063SBarry Smith 
912416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
91317ab2063SBarry Smith {
914416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
915416022c9SBarry Smith   *m = 0; *n = a->m;
91617ab2063SBarry Smith   return 0;
91717ab2063SBarry Smith }
9184e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
91917ab2063SBarry Smith {
920416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
921c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
92217ab2063SBarry Smith 
923416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
92417ab2063SBarry Smith 
925416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
926416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
92717ab2063SBarry Smith   if (idx) {
928416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
9294e093b46SBarry Smith     if (*nz && shift) {
9300452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
93117ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
9324e093b46SBarry Smith     } else if (*nz) {
9334e093b46SBarry Smith       *idx = itmp;
93417ab2063SBarry Smith     }
93517ab2063SBarry Smith     else *idx = 0;
93617ab2063SBarry Smith   }
93717ab2063SBarry Smith   return 0;
93817ab2063SBarry Smith }
93917ab2063SBarry Smith 
9404e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
94117ab2063SBarry Smith {
9424e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
9434e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
94417ab2063SBarry Smith   return 0;
94517ab2063SBarry Smith }
94617ab2063SBarry Smith 
947cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
94817ab2063SBarry Smith {
949416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
950416022c9SBarry Smith   Scalar     *v = a->a;
95117ab2063SBarry Smith   double     sum = 0.0;
952416022c9SBarry Smith   int        i, j,shift = a->indexshift;
95317ab2063SBarry Smith 
95417ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
955416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
95617ab2063SBarry Smith #if defined(PETSC_COMPLEX)
95717ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
95817ab2063SBarry Smith #else
95917ab2063SBarry Smith       sum += (*v)*(*v); v++;
96017ab2063SBarry Smith #endif
96117ab2063SBarry Smith     }
96217ab2063SBarry Smith     *norm = sqrt(sum);
96317ab2063SBarry Smith   }
96417ab2063SBarry Smith   else if (type == NORM_1) {
96517ab2063SBarry Smith     double *tmp;
966416022c9SBarry Smith     int    *jj = a->j;
9670452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
968cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
96917ab2063SBarry Smith     *norm = 0.0;
970416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
971a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
97217ab2063SBarry Smith     }
973416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
97417ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
97517ab2063SBarry Smith     }
9760452661fSBarry Smith     PetscFree(tmp);
97717ab2063SBarry Smith   }
97817ab2063SBarry Smith   else if (type == NORM_INFINITY) {
97917ab2063SBarry Smith     *norm = 0.0;
980416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
981416022c9SBarry Smith       v = a->a + a->i[j] + shift;
98217ab2063SBarry Smith       sum = 0.0;
983416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
984cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
98517ab2063SBarry Smith       }
98617ab2063SBarry Smith       if (sum > *norm) *norm = sum;
98717ab2063SBarry Smith     }
98817ab2063SBarry Smith   }
98917ab2063SBarry Smith   else {
99048d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
99117ab2063SBarry Smith   }
99217ab2063SBarry Smith   return 0;
99317ab2063SBarry Smith }
99417ab2063SBarry Smith 
995416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
99617ab2063SBarry Smith {
997416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
998416022c9SBarry Smith   Mat        C;
999416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1000416022c9SBarry Smith   int        shift = a->indexshift;
1001d5d45c9bSBarry Smith   Scalar     *array = a->a;
100217ab2063SBarry Smith 
10033638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
1004dfa27b74SSatish Balay     SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place");
10050452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1006cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
100717ab2063SBarry Smith   if (shift) {
100817ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
100917ab2063SBarry Smith   }
101017ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1011416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
10120452661fSBarry Smith   PetscFree(col);
101317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
101417ab2063SBarry Smith     len = ai[i+1]-ai[i];
1015416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
101617ab2063SBarry Smith     array += len; aj += len;
101717ab2063SBarry Smith   }
101817ab2063SBarry Smith   if (shift) {
101917ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
102017ab2063SBarry Smith   }
102117ab2063SBarry Smith 
1022416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1023416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
102417ab2063SBarry Smith 
10253638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1026416022c9SBarry Smith     *B = C;
102717ab2063SBarry Smith   } else {
1028416022c9SBarry Smith     /* This isn't really an in-place transpose */
10290452661fSBarry Smith     PetscFree(a->a);
10300452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
10310452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
10320452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
10330452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
10340452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
10351073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
10360452661fSBarry Smith     PetscFree(a);
1037416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
10380452661fSBarry Smith     PetscHeaderDestroy(C);
103917ab2063SBarry Smith   }
104017ab2063SBarry Smith   return 0;
104117ab2063SBarry Smith }
104217ab2063SBarry Smith 
1043f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
104417ab2063SBarry Smith {
1045416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
104617ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1047d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
104817ab2063SBarry Smith 
104917ab2063SBarry Smith   if (ll) {
105017ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
1051f0b747eeSBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
1052416022c9SBarry Smith     v = a->a;
105317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
105417ab2063SBarry Smith       x = l[i];
1055416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
105617ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
105717ab2063SBarry Smith     }
105844cd7ae7SLois Curfman McInnes     PLogFlops(nz);
105917ab2063SBarry Smith   }
106017ab2063SBarry Smith   if (rr) {
106117ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
1062f0b747eeSBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
1063416022c9SBarry Smith     v = a->a; jj = a->j;
106417ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
106517ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
106617ab2063SBarry Smith     }
106744cd7ae7SLois Curfman McInnes     PLogFlops(nz);
106817ab2063SBarry Smith   }
106917ab2063SBarry Smith   return 0;
107017ab2063SBarry Smith }
107117ab2063SBarry Smith 
1072cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
107317ab2063SBarry Smith {
1074db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
107502834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
107699141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1077a2744918SBarry Smith   register int sum,lensi;
107899141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
107999141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
108099141d43SSatish Balay   Scalar       *a_new,*mat_a;
1081416022c9SBarry Smith   Mat          C;
108217ab2063SBarry Smith 
108399141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
108499141d43SSatish Balay   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IS is not sorted");
108599141d43SSatish Balay 
108617ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
108717ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
108817ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
108917ab2063SBarry Smith 
10907264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
109102834360SBarry Smith     /* special case of contiguous rows */
109257faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
109302834360SBarry Smith     starts = lens + ncols;
109402834360SBarry Smith     /* loop over new rows determining lens and starting points */
109502834360SBarry Smith     for (i=0; i<nrows; i++) {
1096a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1097a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
109802834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1099d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
110002834360SBarry Smith           starts[i] = k;
110102834360SBarry Smith           break;
110202834360SBarry Smith 	}
110302834360SBarry Smith       }
1104a2744918SBarry Smith       sum = 0;
110502834360SBarry Smith       while (k < kend) {
1106d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1107a2744918SBarry Smith         sum++;
110802834360SBarry Smith       }
1109a2744918SBarry Smith       lens[i] = sum;
111002834360SBarry Smith     }
111102834360SBarry Smith     /* create submatrix */
1112cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
111308480c60SBarry Smith       int n_cols,n_rows;
111408480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
111508480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1116d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
111708480c60SBarry Smith       C = *B;
111808480c60SBarry Smith     }
111908480c60SBarry Smith     else {
112002834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
112108480c60SBarry Smith     }
1122db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1123db02288aSLois Curfman McInnes 
112402834360SBarry Smith     /* loop over rows inserting into submatrix */
1125db02288aSLois Curfman McInnes     a_new    = c->a;
1126db02288aSLois Curfman McInnes     j_new    = c->j;
1127db02288aSLois Curfman McInnes     i_new    = c->i;
1128db02288aSLois Curfman McInnes     i_new[0] = -shift;
112902834360SBarry Smith     for (i=0; i<nrows; i++) {
1130a2744918SBarry Smith       ii    = starts[i];
1131a2744918SBarry Smith       lensi = lens[i];
1132a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1133a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
113402834360SBarry Smith       }
1135a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1136a2744918SBarry Smith       a_new      += lensi;
1137a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1138a2744918SBarry Smith       c->ilen[i]  = lensi;
113902834360SBarry Smith     }
11400452661fSBarry Smith     PetscFree(lens);
114102834360SBarry Smith   }
114202834360SBarry Smith   else {
114302834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
11440452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
114502834360SBarry Smith     ssmap = smap + shift;
114699141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1147cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
114817ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
114902834360SBarry Smith     /* determine lens of each row */
115002834360SBarry Smith     for (i=0; i<nrows; i++) {
1151d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
115202834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
115302834360SBarry Smith       lens[i] = 0;
115402834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1155d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
115602834360SBarry Smith           lens[i]++;
115702834360SBarry Smith         }
115802834360SBarry Smith       }
115902834360SBarry Smith     }
116017ab2063SBarry Smith     /* Create and fill new matrix */
1161a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
116299141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
116399141d43SSatish Balay 
116499141d43SSatish Balay       if (c->m  != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
116599141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
116699141d43SSatish Balay         SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros");
116799141d43SSatish Balay       }
116899141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
116908480c60SBarry Smith       C = *B;
117008480c60SBarry Smith     }
117108480c60SBarry Smith     else {
117202834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
117308480c60SBarry Smith     }
117499141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
117517ab2063SBarry Smith     for (i=0; i<nrows; i++) {
117699141d43SSatish Balay       row    = irow[i];
117717ab2063SBarry Smith       nznew  = 0;
117899141d43SSatish Balay       kstart = ai[row]+shift;
117999141d43SSatish Balay       kend   = kstart + a->ilen[row];
118099141d43SSatish Balay       mat_i  = c->i[i]+shift;
118199141d43SSatish Balay       mat_j  = c->j + mat_i;
118299141d43SSatish Balay       mat_a  = c->a + mat_i;
118399141d43SSatish Balay       mat_ilen = c->ilen + i;
118417ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
118599141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
118699141d43SSatish Balay           *mat_j++ = tcol - (!shift);
118799141d43SSatish Balay           *mat_a++ = a->a[k];
118899141d43SSatish Balay           (*mat_ilen)++;
118999141d43SSatish Balay 
119017ab2063SBarry Smith         }
119117ab2063SBarry Smith       }
119217ab2063SBarry Smith     }
119302834360SBarry Smith     /* Free work space */
119402834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
119599141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
119602834360SBarry Smith   }
1197416022c9SBarry Smith   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1198416022c9SBarry Smith   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
119917ab2063SBarry Smith 
120017ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1201416022c9SBarry Smith   *B = C;
120217ab2063SBarry Smith   return 0;
120317ab2063SBarry Smith }
120417ab2063SBarry Smith 
1205a871dcd8SBarry Smith /*
120663b91edcSBarry Smith      note: This can only work for identity for row and col. It would
120763b91edcSBarry Smith    be good to check this and otherwise generate an error.
1208a871dcd8SBarry Smith */
120963b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1210a871dcd8SBarry Smith {
121163b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
121208480c60SBarry Smith   int        ierr;
121363b91edcSBarry Smith   Mat        outA;
121463b91edcSBarry Smith 
1215a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1216a871dcd8SBarry Smith 
121763b91edcSBarry Smith   outA          = inA;
121863b91edcSBarry Smith   inA->factor   = FACTOR_LU;
121963b91edcSBarry Smith   a->row        = row;
122063b91edcSBarry Smith   a->col        = col;
122163b91edcSBarry Smith 
12220452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
122363b91edcSBarry Smith 
122408480c60SBarry Smith   if (!a->diag) {
122508480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
122663b91edcSBarry Smith   }
122763b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1228a871dcd8SBarry Smith   return 0;
1229a871dcd8SBarry Smith }
1230a871dcd8SBarry Smith 
1231f0b747eeSBarry Smith #include "pinclude/plapack.h"
1232f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1233f0b747eeSBarry Smith {
1234f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1235f0b747eeSBarry Smith   int        one = 1;
1236f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1237f0b747eeSBarry Smith   PLogFlops(a->nz);
1238f0b747eeSBarry Smith   return 0;
1239f0b747eeSBarry Smith }
1240f0b747eeSBarry Smith 
1241cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1242cddf8d76SBarry Smith                                     Mat **B)
1243cddf8d76SBarry Smith {
1244cddf8d76SBarry Smith   int ierr,i;
1245cddf8d76SBarry Smith 
1246cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
12470452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1248cddf8d76SBarry Smith   }
1249cddf8d76SBarry Smith 
1250cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1251cddf8d76SBarry Smith     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1252cddf8d76SBarry Smith   }
1253cddf8d76SBarry Smith   return 0;
1254cddf8d76SBarry Smith }
1255cddf8d76SBarry Smith 
1256e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
12574dcbc457SBarry Smith {
1258e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
125906763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
12608a047759SSatish Balay   int        start, end, *ai, *aj;
126106763907SSatish Balay   char       *table;
12628a047759SSatish Balay   shift = a->indexshift;
1263e4d965acSSatish Balay   m     = a->m;
1264e4d965acSSatish Balay   ai    = a->i;
12658a047759SSatish Balay   aj    = a->j+shift;
12668a047759SSatish Balay 
1267e4d965acSSatish Balay   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
126806763907SSatish Balay 
126906763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
127006763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
127106763907SSatish Balay 
1272e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1273e4d965acSSatish Balay     /* Initialise the two local arrays */
1274e4d965acSSatish Balay     isz  = 0;
127506763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1276e4d965acSSatish Balay 
1277e4d965acSSatish Balay                 /* Extract the indices, assume there can be duplicate entries */
12784dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
127977c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1280e4d965acSSatish Balay 
1281e4d965acSSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1282e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
128306763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
12844dcbc457SBarry Smith     }
128506763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
128606763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1287e4d965acSSatish Balay 
128804a348a9SBarry Smith     k = 0;
128904a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
129004a348a9SBarry Smith       n = isz;
129106763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1292e4d965acSSatish Balay         row   = nidx[k];
1293e4d965acSSatish Balay         start = ai[row];
1294e4d965acSSatish Balay         end   = ai[row+1];
129504a348a9SBarry Smith         for ( l = start; l<end ; l++){
12968a047759SSatish Balay           val = aj[l] + shift;
129706763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1298e4d965acSSatish Balay         }
1299e4d965acSSatish Balay       }
1300e4d965acSSatish Balay     }
1301e4d965acSSatish Balay     ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1302e4d965acSSatish Balay   }
130304a348a9SBarry Smith   PetscFree(table);
130406763907SSatish Balay   PetscFree(nidx);
1305e4d965acSSatish Balay   return 0;
13064dcbc457SBarry Smith }
130717ab2063SBarry Smith 
1308682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1309682d7d0cSBarry Smith {
1310682d7d0cSBarry Smith   static int called = 0;
1311682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1312682d7d0cSBarry Smith 
1313682d7d0cSBarry Smith   if (called) return 0; else called = 1;
131477c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
131577c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
131677c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
131777c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
131877c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1319682d7d0cSBarry Smith #if defined(HAVE_ESSL)
132077c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1321682d7d0cSBarry Smith #endif
1322682d7d0cSBarry Smith   return 0;
1323682d7d0cSBarry Smith }
1324205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1325682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
132617ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
132717ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1328416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1329416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
133017ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
133117ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
133217ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
133317ab2063SBarry Smith        MatRelax_SeqAIJ,
133417ab2063SBarry Smith        MatTranspose_SeqAIJ,
13357264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1336f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
133717ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
133817ab2063SBarry Smith        MatCompress_SeqAIJ,
133917ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
134017ab2063SBarry Smith        MatGetReordering_SeqAIJ,
134117ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
134217ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
134317ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
134417ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
1345416022c9SBarry Smith        MatGetSubMatrix_SeqAIJ,0,
13463d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1347cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
13487eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1349682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1350f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
1351*a2ce50c7SBarry Smith        MatScale_SeqAIJ,0,0,MatILUDTFactor};
135217ab2063SBarry Smith 
135317ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
135417ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
135517ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
135617ab2063SBarry Smith 
135717ab2063SBarry Smith /*@C
1358682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
13590d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
13606e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
13612bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
13622bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
136317ab2063SBarry Smith 
136417ab2063SBarry Smith    Input Parameters:
136517ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
136617ab2063SBarry Smith .  m - number of rows
136717ab2063SBarry Smith .  n - number of columns
136817ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
13692bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
13702bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
137117ab2063SBarry Smith 
137217ab2063SBarry Smith    Output Parameter:
1373416022c9SBarry Smith .  A - the matrix
137417ab2063SBarry Smith 
137517ab2063SBarry Smith    Notes:
137617ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
137717ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
13780002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
137944cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
138017ab2063SBarry Smith 
138117ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1382a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
13833d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
13843d323bbdSBarry Smith    will get TERRIBLE performance, see the users' manual chapter on
13850d15e28bSLois Curfman McInnes    matrices and the file $(PETSC_DIR)/Performance.
138617ab2063SBarry Smith 
1387682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1388682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1389682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
13906c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
13916c7ebb05SLois Curfman McInnes 
13926c7ebb05SLois Curfman McInnes    Options Database Keys:
13936c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
13946e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
13956e62573dSLois Curfman McInnes $        (max limit=5)
13966e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
13976e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
13986e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
139917ab2063SBarry Smith 
140017ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
140117ab2063SBarry Smith @*/
1402416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
140317ab2063SBarry Smith {
1404416022c9SBarry Smith   Mat        B;
1405416022c9SBarry Smith   Mat_SeqAIJ *b;
140669957df2SSatish Balay   int        i, len, ierr, flg;
1407d5d45c9bSBarry Smith 
1408416022c9SBarry Smith   *A                  = 0;
14090452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1410416022c9SBarry Smith   PLogObjectCreate(B);
14110452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
141244cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1413416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1414416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1415416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1416416022c9SBarry Smith   B->factor           = 0;
1417416022c9SBarry Smith   B->lupivotthreshold = 1.0;
14187a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
141969957df2SSatish Balay                           &flg); CHKERRQ(ierr);
14207a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
14217a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
14227a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1423416022c9SBarry Smith   b->row              = 0;
1424416022c9SBarry Smith   b->col              = 0;
1425416022c9SBarry Smith   b->indexshift       = 0;
1426b810aeb4SBarry Smith   b->reallocs         = 0;
142769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
142869957df2SSatish Balay   if (flg) b->indexshift = -1;
142917ab2063SBarry Smith 
143044cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
143144cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
14320452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1433b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
14347b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
14357b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1436416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
143717ab2063SBarry Smith     nz = nz*m;
143817ab2063SBarry Smith   }
143917ab2063SBarry Smith   else {
144017ab2063SBarry Smith     nz = 0;
1441416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
144217ab2063SBarry Smith   }
144317ab2063SBarry Smith 
144417ab2063SBarry Smith   /* allocate the matrix space */
1445416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
14460452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1447416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1448cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1449416022c9SBarry Smith   b->i  = b->j + nz;
1450416022c9SBarry Smith   b->singlemalloc = 1;
145117ab2063SBarry Smith 
1452416022c9SBarry Smith   b->i[0] = -b->indexshift;
145317ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1454416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
145517ab2063SBarry Smith   }
145617ab2063SBarry Smith 
1457416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
14580452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1459416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1460416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
146117ab2063SBarry Smith 
1462416022c9SBarry Smith   b->nz               = 0;
1463416022c9SBarry Smith   b->maxnz            = nz;
1464416022c9SBarry Smith   b->sorted           = 0;
1465416022c9SBarry Smith   b->roworiented      = 1;
1466416022c9SBarry Smith   b->nonew            = 0;
1467416022c9SBarry Smith   b->diag             = 0;
1468416022c9SBarry Smith   b->solve_work       = 0;
146971bd300dSLois Curfman McInnes   b->spptr            = 0;
1470754ec7b1SSatish Balay   b->inode.node_count = 0;
1471754ec7b1SSatish Balay   b->inode.size       = 0;
14726c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
14736c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
147417ab2063SBarry Smith 
1475416022c9SBarry Smith   *A = B;
14764b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
14774b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
147869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
147969957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
14804b14c69eSBarry Smith #endif
148169957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
148269957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
148369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
148469957df2SSatish Balay   if (flg) {
1485416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1486416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
148717ab2063SBarry Smith   }
148869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
148969957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
149017ab2063SBarry Smith   return 0;
149117ab2063SBarry Smith }
149217ab2063SBarry Smith 
14933d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
149417ab2063SBarry Smith {
1495416022c9SBarry Smith   Mat        C;
1496416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
149708480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
149817ab2063SBarry Smith 
14994043dd9cSLois Curfman McInnes   *B = 0;
15000452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1501416022c9SBarry Smith   PLogObjectCreate(C);
15020452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
150341c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1504416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1505416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1506416022c9SBarry Smith   C->factor     = A->factor;
1507416022c9SBarry Smith   c->row        = 0;
1508416022c9SBarry Smith   c->col        = 0;
1509416022c9SBarry Smith   c->indexshift = shift;
1510c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
151117ab2063SBarry Smith 
151244cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
151344cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
151444cd7ae7SLois Curfman McInnes   C->M          = a->m;
151544cd7ae7SLois Curfman McInnes   C->N          = a->n;
151617ab2063SBarry Smith 
15170452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
15180452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
151917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1520416022c9SBarry Smith     c->imax[i] = a->imax[i];
1521416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
152217ab2063SBarry Smith   }
152317ab2063SBarry Smith 
152417ab2063SBarry Smith   /* allocate the matrix space */
1525416022c9SBarry Smith   c->singlemalloc = 1;
1526416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
15270452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1528416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1529416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1530416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
153117ab2063SBarry Smith   if (m > 0) {
1532416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
153308480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1534416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
153517ab2063SBarry Smith     }
153608480c60SBarry Smith   }
153717ab2063SBarry Smith 
1538416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1539416022c9SBarry Smith   c->sorted      = a->sorted;
1540416022c9SBarry Smith   c->roworiented = a->roworiented;
1541416022c9SBarry Smith   c->nonew       = a->nonew;
15427a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
154317ab2063SBarry Smith 
1544416022c9SBarry Smith   if (a->diag) {
15450452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1546416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
154717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1548416022c9SBarry Smith       c->diag[i] = a->diag[i];
154917ab2063SBarry Smith     }
155017ab2063SBarry Smith   }
1551416022c9SBarry Smith   else c->diag          = 0;
15526c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
15536c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1554754ec7b1SSatish Balay   if (a->inode.size){
1555754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1556754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1557754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1558754ec7b1SSatish Balay   } else {
1559754ec7b1SSatish Balay     c->inode.size       = 0;
1560754ec7b1SSatish Balay     c->inode.node_count = 0;
1561754ec7b1SSatish Balay   }
1562416022c9SBarry Smith   c->nz                 = a->nz;
1563416022c9SBarry Smith   c->maxnz              = a->maxnz;
1564416022c9SBarry Smith   c->solve_work         = 0;
156576dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1566754ec7b1SSatish Balay 
1567416022c9SBarry Smith   *B = C;
156817ab2063SBarry Smith   return 0;
156917ab2063SBarry Smith }
157017ab2063SBarry Smith 
157119bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
157217ab2063SBarry Smith {
1573416022c9SBarry Smith   Mat_SeqAIJ   *a;
1574416022c9SBarry Smith   Mat          B;
157517699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1576bcd2baecSBarry Smith   MPI_Comm     comm;
157717ab2063SBarry Smith 
157819bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
157917699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
158017699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
158190ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
158277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
158348d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
158417ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
158517ab2063SBarry Smith 
158617ab2063SBarry Smith   /* read in row lengths */
15870452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
158877c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
158917ab2063SBarry Smith 
159017ab2063SBarry Smith   /* create our matrix */
1591416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1592416022c9SBarry Smith   B = *A;
1593416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1594416022c9SBarry Smith   shift = a->indexshift;
159517ab2063SBarry Smith 
159617ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
159777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
159817ab2063SBarry Smith   if (shift) {
159917ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1600416022c9SBarry Smith       a->j[i] += 1;
160117ab2063SBarry Smith     }
160217ab2063SBarry Smith   }
160317ab2063SBarry Smith 
160417ab2063SBarry Smith   /* read in nonzero values */
160577c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
160617ab2063SBarry Smith 
160717ab2063SBarry Smith   /* set matrix "i" values */
1608416022c9SBarry Smith   a->i[0] = -shift;
160917ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1610416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1611416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
161217ab2063SBarry Smith   }
16130452661fSBarry Smith   PetscFree(rowlengths);
161417ab2063SBarry Smith 
1615416022c9SBarry Smith   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1616416022c9SBarry Smith   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
161717ab2063SBarry Smith   return 0;
161817ab2063SBarry Smith }
161917ab2063SBarry Smith 
1620686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
16217264ac53SSatish Balay {
16227264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
16237264ac53SSatish Balay 
1624bcd2baecSBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type");
16257264ac53SSatish Balay 
16267264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
16277264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1628bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
162977c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1630bcd2baecSBarry Smith   }
16317264ac53SSatish Balay 
16327264ac53SSatish Balay   /* if the a->i are the same */
1633bcd2baecSBarry Smith   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
163477c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
16357264ac53SSatish Balay   }
16367264ac53SSatish Balay 
16377264ac53SSatish Balay   /* if a->j are the same */
1638bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
163977c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1640bcd2baecSBarry Smith   }
1641bcd2baecSBarry Smith 
1642bcd2baecSBarry Smith   /* if a->a are the same */
164319bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
164477c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
16457264ac53SSatish Balay   }
164677c4ece6SBarry Smith   *flg = PETSC_TRUE;
16477264ac53SSatish Balay   return 0;
16487264ac53SSatish Balay 
16497264ac53SSatish Balay }
1650