xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 4e220ebcc634eb0c8ccd78bed049803c9194c4d4)
16d84be18SBarry Smith 
217ab2063SBarry Smith #ifndef lint
3*4e220ebcSLois Curfman McInnes static char vcid[] = "$Id: aij.c,v 1.182 1996/08/15 12:47:21 bsmith Exp curfman $";
417ab2063SBarry Smith #endif
517ab2063SBarry Smith 
6d5d45c9bSBarry Smith /*
75a838052SSatish Balay B    Defines the basic matrix operations for the AIJ (compressed row)
8d5d45c9bSBarry Smith   matrix storage format.
9d5d45c9bSBarry Smith */
1070f55243SBarry Smith #include "src/mat/impls/aij/seq/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 
16a2ce50c7SBarry Smith /*
17a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
18a2ce50c7SBarry Smith */
19a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
20a2ce50c7SBarry Smith {
21a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
22a2ce50c7SBarry Smith   int        ierr = 1;
23a2ce50c7SBarry Smith 
24a2ce50c7SBarry Smith   SETERRQ(ierr,"MatILUDTFactor_SeqAIJ:Not implemented");
25a2ce50c7SBarry Smith }
26a2ce50c7SBarry Smith 
27a2ce50c7SBarry Smith /*
28a2ce50c7SBarry Smith      Basic flow based reordering routine for AIJ storage format
29a2ce50c7SBarry Smith */
30a2ce50c7SBarry Smith int MatGetReordering_SeqAIJ_Flow(Mat A,IS *rperm,IS *cperm)
31a2ce50c7SBarry Smith {
32a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
33a2ce50c7SBarry Smith   int        ierr = 1;
34a2ce50c7SBarry Smith 
35a2ce50c7SBarry Smith   SETERRQ(ierr,"MatGetReordering_SeqAIJ_Flow:Not implemented");
36a2ce50c7SBarry Smith }
37a2ce50c7SBarry Smith 
38bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3917ab2063SBarry Smith 
40a2ce50c7SBarry 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   /*
46a2ce50c7SBarry Smith      This is tacky. The problem is we cannot use the registered ordering
47a2ce50c7SBarry Smith     routines since they only work on nozero patterns. To have general
48a2ce50c7SBarry Smith     registration means you register a particular ordering method for a
49a2ce50c7SBarry Smith     particular data structure; hence something like a two dimensional
50a2ce50c7SBarry Smith     look up table.
51a2ce50c7SBarry Smith   */
52a2ce50c7SBarry Smith   if (type == ORDER_FLOW) {
53a2ce50c7SBarry Smith      return MatGetReordering_SeqAIJ_Flow(A,rperm,cperm);
54a2ce50c7SBarry Smith   }
55a2ce50c7SBarry Smith 
56a2ce50c7SBarry 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;
65537820f0SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
66537820f0SBarry Smith     ierr = ISCreateGeneral(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) {
295416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
296*4e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
297*4e220ebcSLois Curfman McInnes     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz);
29817ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
29917ab2063SBarry Smith 
30017ab2063SBarry Smith     for (i=0; i<m; i++) {
301416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
30217ab2063SBarry Smith #if defined(PETSC_COMPLEX)
3037a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]),
304416022c9SBarry Smith                    imag(a->a[j]));
30517ab2063SBarry Smith #else
3067a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
30717ab2063SBarry Smith #endif
30817ab2063SBarry Smith       }
30917ab2063SBarry Smith     }
31017ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
31117ab2063SBarry Smith   }
31244cd7ae7SLois Curfman McInnes   else if (format == ASCII_FORMAT_COMMON) {
31344cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
31444cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
31544cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
31644cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
31744cd7ae7SLois Curfman McInnes         if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0)
31844cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
31944cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
32044cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
32144cd7ae7SLois Curfman McInnes #else
32244cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
32344cd7ae7SLois Curfman McInnes #endif
32444cd7ae7SLois Curfman McInnes       }
32544cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
32644cd7ae7SLois Curfman McInnes     }
32744cd7ae7SLois Curfman McInnes   }
32817ab2063SBarry Smith   else {
32917ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
33017ab2063SBarry Smith       fprintf(fd,"row %d:",i);
331416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
33217ab2063SBarry Smith #if defined(PETSC_COMPLEX)
333416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
334416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
33517ab2063SBarry Smith         }
33617ab2063SBarry Smith         else {
337416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
33817ab2063SBarry Smith         }
33917ab2063SBarry Smith #else
340416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
34117ab2063SBarry Smith #endif
34217ab2063SBarry Smith       }
34317ab2063SBarry Smith       fprintf(fd,"\n");
34417ab2063SBarry Smith     }
34517ab2063SBarry Smith   }
34617ab2063SBarry Smith   fflush(fd);
347416022c9SBarry Smith   return 0;
348416022c9SBarry Smith }
349416022c9SBarry Smith 
350416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
351416022c9SBarry Smith {
352416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
353cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
354cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
355bcd2baecSBarry Smith   Draw        draw;
356cddf8d76SBarry Smith   DrawButton  button;
35719bcc07fSBarry Smith   PetscTruth  isnull;
358cddf8d76SBarry Smith 
359bcd2baecSBarry Smith   ViewerDrawGetDraw(viewer,&draw);
36019bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
36119bcc07fSBarry Smith 
362416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
363416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
364416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
365416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
366cddf8d76SBarry Smith   color = DRAW_BLUE;
367416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
368cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
369416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
370cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
371cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
372cddf8d76SBarry Smith       if (real(a->a[j]) >=  0.) continue;
373cddf8d76SBarry Smith #else
374cddf8d76SBarry Smith       if (a->a[j] >=  0.) continue;
375cddf8d76SBarry Smith #endif
376cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
377cddf8d76SBarry Smith     }
378cddf8d76SBarry Smith   }
379cddf8d76SBarry Smith   color = DRAW_CYAN;
380cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
381cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
382cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
383cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
384cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
385cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
386cddf8d76SBarry Smith     }
387cddf8d76SBarry Smith   }
388cddf8d76SBarry Smith   color = DRAW_RED;
389cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
390cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
391cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
392cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
393cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
394cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
395cddf8d76SBarry Smith #else
396cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
397cddf8d76SBarry Smith #endif
398cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
399416022c9SBarry Smith     }
400416022c9SBarry Smith   }
401416022c9SBarry Smith   DrawFlush(draw);
402cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
403cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
404cddf8d76SBarry Smith 
405cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
406cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
407cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
408cddf8d76SBarry Smith     DrawClear(draw);
409cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
410cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
411cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
412cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
413cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
414cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
415cddf8d76SBarry Smith     w *= scale; h *= scale;
416cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
417cddf8d76SBarry Smith     color = DRAW_BLUE;
418cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
419cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
420cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
421cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
422cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
423cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
424cddf8d76SBarry Smith #else
425cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
426cddf8d76SBarry Smith #endif
427cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
428cddf8d76SBarry Smith       }
429cddf8d76SBarry Smith     }
430cddf8d76SBarry Smith     color = DRAW_CYAN;
431cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
432cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
433cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
434cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
435cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
436cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
437cddf8d76SBarry Smith       }
438cddf8d76SBarry Smith     }
439cddf8d76SBarry Smith     color = DRAW_RED;
440cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
441cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
442cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
443cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
444cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
445cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
446cddf8d76SBarry Smith #else
447cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
448cddf8d76SBarry Smith #endif
449cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
450cddf8d76SBarry Smith       }
451cddf8d76SBarry Smith     }
452cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
453cddf8d76SBarry Smith   }
454416022c9SBarry Smith   return 0;
455416022c9SBarry Smith }
456416022c9SBarry Smith 
457416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
458416022c9SBarry Smith {
459416022c9SBarry Smith   Mat         A = (Mat) obj;
460416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
461bcd2baecSBarry Smith   ViewerType  vtype;
462bcd2baecSBarry Smith   int         ierr;
463416022c9SBarry Smith 
464bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
465bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
466416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
467416022c9SBarry Smith   }
468bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
469416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
470416022c9SBarry Smith   }
471bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
472416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
473416022c9SBarry Smith   }
474bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
475bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
47617ab2063SBarry Smith   }
47717ab2063SBarry Smith   return 0;
47817ab2063SBarry Smith }
47919bcc07fSBarry Smith 
480c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
481416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
48217ab2063SBarry Smith {
483416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
48441c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
485416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
486416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
48717ab2063SBarry Smith 
4886d4a8577SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
48917ab2063SBarry Smith 
49017ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
491416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
49217ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
49317ab2063SBarry Smith     if (fshift) {
494416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
49517ab2063SBarry Smith       N = ailen[i];
49617ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
49717ab2063SBarry Smith         ip[j-fshift] = ip[j];
49817ab2063SBarry Smith         ap[j-fshift] = ap[j];
49917ab2063SBarry Smith       }
50017ab2063SBarry Smith     }
50117ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
50217ab2063SBarry Smith   }
50317ab2063SBarry Smith   if (m) {
50417ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
50517ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
50617ab2063SBarry Smith   }
50717ab2063SBarry Smith   /* reset ilen and imax for each row */
50817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
50917ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
51017ab2063SBarry Smith   }
511416022c9SBarry Smith   a->nz = ai[m] + shift;
51217ab2063SBarry Smith 
51317ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
514416022c9SBarry Smith   if (fshift && a->diag) {
5150452661fSBarry Smith     PetscFree(a->diag);
516416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
517416022c9SBarry Smith     a->diag = 0;
51817ab2063SBarry Smith   }
519*4e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
520*4e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
521*4e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
522b810aeb4SBarry Smith            a->reallocs);
523*4e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
524*4e220ebcSLois Curfman McInnes 
52576dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
52641c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
52717ab2063SBarry Smith   return 0;
52817ab2063SBarry Smith }
52917ab2063SBarry Smith 
530416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
53117ab2063SBarry Smith {
532416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
533cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
53417ab2063SBarry Smith   return 0;
53517ab2063SBarry Smith }
536416022c9SBarry Smith 
53717ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
53817ab2063SBarry Smith {
539416022c9SBarry Smith   Mat        A  = (Mat) obj;
540416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
541d5d45c9bSBarry Smith 
54217ab2063SBarry Smith #if defined(PETSC_LOG)
543416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
54417ab2063SBarry Smith #endif
5450452661fSBarry Smith   PetscFree(a->a);
5460452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
5470452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
5480452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
5490452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
5500452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
55176dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
5520452661fSBarry Smith   PetscFree(a);
553f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
554f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
55517ab2063SBarry Smith   return 0;
55617ab2063SBarry Smith }
55717ab2063SBarry Smith 
558416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
55917ab2063SBarry Smith {
56017ab2063SBarry Smith   return 0;
56117ab2063SBarry Smith }
56217ab2063SBarry Smith 
563416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
56417ab2063SBarry Smith {
565416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
5666d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)              a->roworiented = 1;
5676d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)           a->roworiented = 0;
5686d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)            a->sorted      = 1;
5696d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
5706d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
5716d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
5726d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
5736d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
5746d4a8577SBarry Smith            op == MAT_YES_NEW_DIAGONALS)
57594a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
5766d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
5776d4a8577SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:MAT_NO_NEW_DIAGONALS");}
5786d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
5796d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
5806d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
5816d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
5826d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
583e2f28af5SBarry Smith   else
5844d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
58517ab2063SBarry Smith   return 0;
58617ab2063SBarry Smith }
58717ab2063SBarry Smith 
588416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
58917ab2063SBarry Smith {
590416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
591416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
59217ab2063SBarry Smith   Scalar     *x, zero = 0.0;
59317ab2063SBarry Smith 
59417ab2063SBarry Smith   VecSet(&zero,v);
59517ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
596416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
597416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
598416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
599416022c9SBarry Smith       if (a->j[j]+shift == i) {
600416022c9SBarry Smith         x[i] = a->a[j];
60117ab2063SBarry Smith         break;
60217ab2063SBarry Smith       }
60317ab2063SBarry Smith     }
60417ab2063SBarry Smith   }
60517ab2063SBarry Smith   return 0;
60617ab2063SBarry Smith }
60717ab2063SBarry Smith 
60817ab2063SBarry Smith /* -------------------------------------------------------*/
60917ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
61017ab2063SBarry Smith /* -------------------------------------------------------*/
61144cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
61217ab2063SBarry Smith {
613416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
61417ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
615416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
61617ab2063SBarry Smith 
61717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
618cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
61917ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
62017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
621416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
622416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
623416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
62417ab2063SBarry Smith     alpha = x[i];
62517ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
62617ab2063SBarry Smith   }
627416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
62817ab2063SBarry Smith   return 0;
62917ab2063SBarry Smith }
630d5d45c9bSBarry Smith 
63144cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
63217ab2063SBarry Smith {
633416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
63417ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
635416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
63617ab2063SBarry Smith 
63717ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
63817ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
63917ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
64017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
641416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
642416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
643416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
64417ab2063SBarry Smith     alpha = x[i];
64517ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
64617ab2063SBarry Smith   }
64717ab2063SBarry Smith   return 0;
64817ab2063SBarry Smith }
64917ab2063SBarry Smith 
65044cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
65117ab2063SBarry Smith {
652416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
65317ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
654416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
65517ab2063SBarry Smith 
65617ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
65717ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
658416022c9SBarry Smith   idx  = a->j;
659416022c9SBarry Smith   v    = a->a;
660416022c9SBarry Smith   ii   = a->i;
66117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
662416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
66317ab2063SBarry Smith     sum  = 0.0;
66417ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
66517ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
666416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
66717ab2063SBarry Smith     y[i] = sum;
66817ab2063SBarry Smith   }
669416022c9SBarry Smith   PLogFlops(2*a->nz - m);
67017ab2063SBarry Smith   return 0;
67117ab2063SBarry Smith }
67217ab2063SBarry Smith 
67344cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
67417ab2063SBarry Smith {
675416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
67617ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
677cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
67817ab2063SBarry Smith 
67917ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
68017ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
681cddf8d76SBarry Smith   idx  = a->j;
682cddf8d76SBarry Smith   v    = a->a;
683cddf8d76SBarry Smith   ii   = a->i;
68417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
685cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
68617ab2063SBarry Smith     sum  = y[i];
687cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
688cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
68917ab2063SBarry Smith     z[i] = sum;
69017ab2063SBarry Smith   }
691416022c9SBarry Smith   PLogFlops(2*a->nz);
69217ab2063SBarry Smith   return 0;
69317ab2063SBarry Smith }
69417ab2063SBarry Smith 
69517ab2063SBarry Smith /*
69617ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
69717ab2063SBarry Smith */
69817ab2063SBarry Smith 
699416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
70017ab2063SBarry Smith {
701416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
702416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
70317ab2063SBarry Smith 
7040452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
705416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
706416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
707416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
708416022c9SBarry Smith       if (a->j[j]+shift == i) {
70917ab2063SBarry Smith         diag[i] = j - shift;
71017ab2063SBarry Smith         break;
71117ab2063SBarry Smith       }
71217ab2063SBarry Smith     }
71317ab2063SBarry Smith   }
714416022c9SBarry Smith   a->diag = diag;
71517ab2063SBarry Smith   return 0;
71617ab2063SBarry Smith }
71717ab2063SBarry Smith 
71844cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
71917ab2063SBarry Smith                            double fshift,int its,Vec xx)
72017ab2063SBarry Smith {
721416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
722416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
723d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
72417ab2063SBarry Smith 
72517ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
726416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
727416022c9SBarry Smith   diag = a->diag;
728416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
72917ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
73017ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
73117ab2063SBarry Smith     bs = b + shift;
73217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
733416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
734416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
735416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
736416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
73717ab2063SBarry Smith         sum  = b[i]*d/omega;
73817ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
73917ab2063SBarry Smith         x[i] = sum;
74017ab2063SBarry Smith     }
74117ab2063SBarry Smith     return 0;
74217ab2063SBarry Smith   }
74317ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
74417ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
74517ab2063SBarry Smith   }
746416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
74717ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
74817ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
74917ab2063SBarry Smith 
75017ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
75117ab2063SBarry Smith 
75217ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
75317ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
75417ab2063SBarry Smith     is the relaxation factor.
75517ab2063SBarry Smith     */
7560452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
75717ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
75817ab2063SBarry Smith 
75917ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
76017ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
761416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
762416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
763416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
764416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
76517ab2063SBarry Smith       sum  = b[i];
76617ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
76717ab2063SBarry Smith       x[i] = omega*(sum/d);
76817ab2063SBarry Smith     }
76917ab2063SBarry Smith 
77017ab2063SBarry Smith     /*  t = b - (2*E - D)x */
771416022c9SBarry Smith     v = a->a;
77217ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
77317ab2063SBarry Smith 
77417ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
775416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
776416022c9SBarry Smith     diag = a->diag;
77717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
778416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
779416022c9SBarry Smith       n    = diag[i] - a->i[i];
780416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
781416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
78217ab2063SBarry Smith       sum  = t[i];
78317ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
78417ab2063SBarry Smith       t[i] = omega*(sum/d);
78517ab2063SBarry Smith     }
78617ab2063SBarry Smith 
78717ab2063SBarry Smith     /*  x = x + t */
78817ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
7890452661fSBarry Smith     PetscFree(t);
79017ab2063SBarry Smith     return 0;
79117ab2063SBarry Smith   }
79217ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
79317ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
79417ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
795416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
796416022c9SBarry Smith         n    = diag[i] - a->i[i];
797416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
798416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
79917ab2063SBarry Smith         sum  = b[i];
80017ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
80117ab2063SBarry Smith         x[i] = omega*(sum/d);
80217ab2063SBarry Smith       }
80317ab2063SBarry Smith       xb = x;
80417ab2063SBarry Smith     }
80517ab2063SBarry Smith     else xb = b;
80617ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
80717ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
80817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
809416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
81017ab2063SBarry Smith       }
81117ab2063SBarry Smith     }
81217ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
81317ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
814416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
815416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
816416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
817416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
81817ab2063SBarry Smith         sum  = xb[i];
81917ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
82017ab2063SBarry Smith         x[i] = omega*(sum/d);
82117ab2063SBarry Smith       }
82217ab2063SBarry Smith     }
82317ab2063SBarry Smith     its--;
82417ab2063SBarry Smith   }
82517ab2063SBarry Smith   while (its--) {
82617ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
82717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
828416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
829416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
830416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
831416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
83217ab2063SBarry Smith         sum  = b[i];
83317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
83417ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
83517ab2063SBarry Smith       }
83617ab2063SBarry Smith     }
83717ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
83817ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
839416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
840416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
841416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
842416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
84317ab2063SBarry Smith         sum  = b[i];
84417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
84517ab2063SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
84617ab2063SBarry Smith       }
84717ab2063SBarry Smith     }
84817ab2063SBarry Smith   }
84917ab2063SBarry Smith   return 0;
85017ab2063SBarry Smith }
85117ab2063SBarry Smith 
852*4e220ebcSLois Curfman McInnes static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
85317ab2063SBarry Smith {
854416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
855*4e220ebcSLois Curfman McInnes 
856*4e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
857*4e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
858*4e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
859*4e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
860*4e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
861*4e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
862*4e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
863*4e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
864*4e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
865*4e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
866*4e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
867*4e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
868*4e220ebcSLois Curfman McInnes   info->memory         = A->mem;
869*4e220ebcSLois Curfman McInnes   if (A->factor) {
870*4e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
871*4e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
872*4e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
873*4e220ebcSLois Curfman McInnes   } else {
874*4e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
875*4e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
876*4e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
877*4e220ebcSLois Curfman McInnes   }
87817ab2063SBarry Smith   return 0;
87917ab2063SBarry Smith }
88017ab2063SBarry Smith 
88117ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
88217ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
88317ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
88417ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
88517ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
88617ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
88717ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
88817ab2063SBarry Smith 
88917ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
89017ab2063SBarry Smith {
891416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
892416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
89317ab2063SBarry Smith 
89477c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
89517ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
89617ab2063SBarry Smith   if (diag) {
89717ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
898416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
899416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
900416022c9SBarry Smith         a->ilen[rows[i]] = 1;
901416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
902416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
90317ab2063SBarry Smith       }
90417ab2063SBarry Smith       else {
90517ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
90617ab2063SBarry Smith         CHKERRQ(ierr);
90717ab2063SBarry Smith       }
90817ab2063SBarry Smith     }
90917ab2063SBarry Smith   }
91017ab2063SBarry Smith   else {
91117ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
912416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
913416022c9SBarry Smith       a->ilen[rows[i]] = 0;
91417ab2063SBarry Smith     }
91517ab2063SBarry Smith   }
91617ab2063SBarry Smith   ISRestoreIndices(is,&rows);
9176d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9186d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
91917ab2063SBarry Smith   return 0;
92017ab2063SBarry Smith }
92117ab2063SBarry Smith 
922416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
92317ab2063SBarry Smith {
924416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
925416022c9SBarry Smith   *m = a->m; *n = a->n;
92617ab2063SBarry Smith   return 0;
92717ab2063SBarry Smith }
92817ab2063SBarry Smith 
929416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
93017ab2063SBarry Smith {
931416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
932416022c9SBarry Smith   *m = 0; *n = a->m;
93317ab2063SBarry Smith   return 0;
93417ab2063SBarry Smith }
9354e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
93617ab2063SBarry Smith {
937416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
938c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
93917ab2063SBarry Smith 
940416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
94117ab2063SBarry Smith 
942416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
943416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
94417ab2063SBarry Smith   if (idx) {
945416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
9464e093b46SBarry Smith     if (*nz && shift) {
9470452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
94817ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
9494e093b46SBarry Smith     } else if (*nz) {
9504e093b46SBarry Smith       *idx = itmp;
95117ab2063SBarry Smith     }
95217ab2063SBarry Smith     else *idx = 0;
95317ab2063SBarry Smith   }
95417ab2063SBarry Smith   return 0;
95517ab2063SBarry Smith }
95617ab2063SBarry Smith 
9574e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
95817ab2063SBarry Smith {
9594e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
9604e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
96117ab2063SBarry Smith   return 0;
96217ab2063SBarry Smith }
96317ab2063SBarry Smith 
964cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
96517ab2063SBarry Smith {
966416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
967416022c9SBarry Smith   Scalar     *v = a->a;
96817ab2063SBarry Smith   double     sum = 0.0;
969416022c9SBarry Smith   int        i, j,shift = a->indexshift;
97017ab2063SBarry Smith 
97117ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
972416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
97317ab2063SBarry Smith #if defined(PETSC_COMPLEX)
97417ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
97517ab2063SBarry Smith #else
97617ab2063SBarry Smith       sum += (*v)*(*v); v++;
97717ab2063SBarry Smith #endif
97817ab2063SBarry Smith     }
97917ab2063SBarry Smith     *norm = sqrt(sum);
98017ab2063SBarry Smith   }
98117ab2063SBarry Smith   else if (type == NORM_1) {
98217ab2063SBarry Smith     double *tmp;
983416022c9SBarry Smith     int    *jj = a->j;
9840452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
985cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
98617ab2063SBarry Smith     *norm = 0.0;
987416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
988a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
98917ab2063SBarry Smith     }
990416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
99117ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
99217ab2063SBarry Smith     }
9930452661fSBarry Smith     PetscFree(tmp);
99417ab2063SBarry Smith   }
99517ab2063SBarry Smith   else if (type == NORM_INFINITY) {
99617ab2063SBarry Smith     *norm = 0.0;
997416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
998416022c9SBarry Smith       v = a->a + a->i[j] + shift;
99917ab2063SBarry Smith       sum = 0.0;
1000416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1001cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
100217ab2063SBarry Smith       }
100317ab2063SBarry Smith       if (sum > *norm) *norm = sum;
100417ab2063SBarry Smith     }
100517ab2063SBarry Smith   }
100617ab2063SBarry Smith   else {
100748d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
100817ab2063SBarry Smith   }
100917ab2063SBarry Smith   return 0;
101017ab2063SBarry Smith }
101117ab2063SBarry Smith 
1012416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
101317ab2063SBarry Smith {
1014416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1015416022c9SBarry Smith   Mat        C;
1016416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1017416022c9SBarry Smith   int        shift = a->indexshift;
1018d5d45c9bSBarry Smith   Scalar     *array = a->a;
101917ab2063SBarry Smith 
10203638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
1021dfa27b74SSatish Balay     SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place");
10220452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1023cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
102417ab2063SBarry Smith   if (shift) {
102517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
102617ab2063SBarry Smith   }
102717ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1028416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
10290452661fSBarry Smith   PetscFree(col);
103017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
103117ab2063SBarry Smith     len = ai[i+1]-ai[i];
1032416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
103317ab2063SBarry Smith     array += len; aj += len;
103417ab2063SBarry Smith   }
103517ab2063SBarry Smith   if (shift) {
103617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
103717ab2063SBarry Smith   }
103817ab2063SBarry Smith 
10396d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10406d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
104117ab2063SBarry Smith 
10423638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1043416022c9SBarry Smith     *B = C;
104417ab2063SBarry Smith   } else {
1045416022c9SBarry Smith     /* This isn't really an in-place transpose */
10460452661fSBarry Smith     PetscFree(a->a);
10470452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
10480452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
10490452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
10500452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
10510452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
10521073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
10530452661fSBarry Smith     PetscFree(a);
1054416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
10550452661fSBarry Smith     PetscHeaderDestroy(C);
105617ab2063SBarry Smith   }
105717ab2063SBarry Smith   return 0;
105817ab2063SBarry Smith }
105917ab2063SBarry Smith 
1060f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
106117ab2063SBarry Smith {
1062416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
106317ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1064d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
106517ab2063SBarry Smith 
106617ab2063SBarry Smith   if (ll) {
106717ab2063SBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
1068f0b747eeSBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
1069416022c9SBarry Smith     v = a->a;
107017ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
107117ab2063SBarry Smith       x = l[i];
1072416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
107317ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
107417ab2063SBarry Smith     }
107544cd7ae7SLois Curfman McInnes     PLogFlops(nz);
107617ab2063SBarry Smith   }
107717ab2063SBarry Smith   if (rr) {
107817ab2063SBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
1079f0b747eeSBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
1080416022c9SBarry Smith     v = a->a; jj = a->j;
108117ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
108217ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
108317ab2063SBarry Smith     }
108444cd7ae7SLois Curfman McInnes     PLogFlops(nz);
108517ab2063SBarry Smith   }
108617ab2063SBarry Smith   return 0;
108717ab2063SBarry Smith }
108817ab2063SBarry Smith 
1089cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
109017ab2063SBarry Smith {
1091db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
109202834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
109399141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1094a2744918SBarry Smith   register int sum,lensi;
109599141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
109699141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
109799141d43SSatish Balay   Scalar       *a_new,*mat_a;
1098416022c9SBarry Smith   Mat          C;
109917ab2063SBarry Smith 
1100b48a1e75SSatish Balay   ierr = ISSorted(isrow,(PetscTruth*)&i);
1101b48a1e75SSatish Balay   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:ISrow is not sorted");
110299141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
1103b48a1e75SSatish Balay   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IScol is not sorted");
110499141d43SSatish Balay 
110517ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
110617ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
110717ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
110817ab2063SBarry Smith 
11097264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
111002834360SBarry Smith     /* special case of contiguous rows */
111157faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
111202834360SBarry Smith     starts = lens + ncols;
111302834360SBarry Smith     /* loop over new rows determining lens and starting points */
111402834360SBarry Smith     for (i=0; i<nrows; i++) {
1115a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1116a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
111702834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1118d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
111902834360SBarry Smith           starts[i] = k;
112002834360SBarry Smith           break;
112102834360SBarry Smith 	}
112202834360SBarry Smith       }
1123a2744918SBarry Smith       sum = 0;
112402834360SBarry Smith       while (k < kend) {
1125d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1126a2744918SBarry Smith         sum++;
112702834360SBarry Smith       }
1128a2744918SBarry Smith       lens[i] = sum;
112902834360SBarry Smith     }
113002834360SBarry Smith     /* create submatrix */
1131cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
113208480c60SBarry Smith       int n_cols,n_rows;
113308480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
113408480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1135d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
113608480c60SBarry Smith       C = *B;
113708480c60SBarry Smith     }
113808480c60SBarry Smith     else {
113902834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
114008480c60SBarry Smith     }
1141db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1142db02288aSLois Curfman McInnes 
114302834360SBarry Smith     /* loop over rows inserting into submatrix */
1144db02288aSLois Curfman McInnes     a_new    = c->a;
1145db02288aSLois Curfman McInnes     j_new    = c->j;
1146db02288aSLois Curfman McInnes     i_new    = c->i;
1147db02288aSLois Curfman McInnes     i_new[0] = -shift;
114802834360SBarry Smith     for (i=0; i<nrows; i++) {
1149a2744918SBarry Smith       ii    = starts[i];
1150a2744918SBarry Smith       lensi = lens[i];
1151a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1152a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
115302834360SBarry Smith       }
1154a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1155a2744918SBarry Smith       a_new      += lensi;
1156a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1157a2744918SBarry Smith       c->ilen[i]  = lensi;
115802834360SBarry Smith     }
11590452661fSBarry Smith     PetscFree(lens);
116002834360SBarry Smith   }
116102834360SBarry Smith   else {
116202834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
11630452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
116402834360SBarry Smith     ssmap = smap + shift;
116599141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1166cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
116717ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
116802834360SBarry Smith     /* determine lens of each row */
116902834360SBarry Smith     for (i=0; i<nrows; i++) {
1170d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
117102834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
117202834360SBarry Smith       lens[i] = 0;
117302834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1174d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
117502834360SBarry Smith           lens[i]++;
117602834360SBarry Smith         }
117702834360SBarry Smith       }
117802834360SBarry Smith     }
117917ab2063SBarry Smith     /* Create and fill new matrix */
1180a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
118199141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
118299141d43SSatish Balay 
118399141d43SSatish Balay       if (c->m  != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
118499141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
118599141d43SSatish Balay         SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros");
118699141d43SSatish Balay       }
118799141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
118808480c60SBarry Smith       C = *B;
118908480c60SBarry Smith     }
119008480c60SBarry Smith     else {
119102834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
119208480c60SBarry Smith     }
119399141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
119417ab2063SBarry Smith     for (i=0; i<nrows; i++) {
119599141d43SSatish Balay       row    = irow[i];
119617ab2063SBarry Smith       nznew  = 0;
119799141d43SSatish Balay       kstart = ai[row]+shift;
119899141d43SSatish Balay       kend   = kstart + a->ilen[row];
119999141d43SSatish Balay       mat_i  = c->i[i]+shift;
120099141d43SSatish Balay       mat_j  = c->j + mat_i;
120199141d43SSatish Balay       mat_a  = c->a + mat_i;
120299141d43SSatish Balay       mat_ilen = c->ilen + i;
120317ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
120499141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
120599141d43SSatish Balay           *mat_j++ = tcol - (!shift);
120699141d43SSatish Balay           *mat_a++ = a->a[k];
120799141d43SSatish Balay           (*mat_ilen)++;
120899141d43SSatish Balay 
120917ab2063SBarry Smith         }
121017ab2063SBarry Smith       }
121117ab2063SBarry Smith     }
121202834360SBarry Smith     /* Free work space */
121302834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
121499141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
121502834360SBarry Smith   }
12166d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12176d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
121817ab2063SBarry Smith 
121917ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1220416022c9SBarry Smith   *B = C;
122117ab2063SBarry Smith   return 0;
122217ab2063SBarry Smith }
122317ab2063SBarry Smith 
1224a871dcd8SBarry Smith /*
122563b91edcSBarry Smith      note: This can only work for identity for row and col. It would
122663b91edcSBarry Smith    be good to check this and otherwise generate an error.
1227a871dcd8SBarry Smith */
122863b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1229a871dcd8SBarry Smith {
123063b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
123108480c60SBarry Smith   int        ierr;
123263b91edcSBarry Smith   Mat        outA;
123363b91edcSBarry Smith 
1234a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1235a871dcd8SBarry Smith 
123663b91edcSBarry Smith   outA          = inA;
123763b91edcSBarry Smith   inA->factor   = FACTOR_LU;
123863b91edcSBarry Smith   a->row        = row;
123963b91edcSBarry Smith   a->col        = col;
124063b91edcSBarry Smith 
12410452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
124263b91edcSBarry Smith 
124308480c60SBarry Smith   if (!a->diag) {
124408480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
124563b91edcSBarry Smith   }
124663b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1247a871dcd8SBarry Smith   return 0;
1248a871dcd8SBarry Smith }
1249a871dcd8SBarry Smith 
1250f0b747eeSBarry Smith #include "pinclude/plapack.h"
1251f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1252f0b747eeSBarry Smith {
1253f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1254f0b747eeSBarry Smith   int        one = 1;
1255f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1256f0b747eeSBarry Smith   PLogFlops(a->nz);
1257f0b747eeSBarry Smith   return 0;
1258f0b747eeSBarry Smith }
1259f0b747eeSBarry Smith 
1260cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1261cddf8d76SBarry Smith                                     Mat **B)
1262cddf8d76SBarry Smith {
1263cddf8d76SBarry Smith   int ierr,i;
1264cddf8d76SBarry Smith 
1265cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
12660452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1267cddf8d76SBarry Smith   }
1268cddf8d76SBarry Smith 
1269cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1270905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1271cddf8d76SBarry Smith   }
1272cddf8d76SBarry Smith   return 0;
1273cddf8d76SBarry Smith }
1274cddf8d76SBarry Smith 
12755a838052SSatish Balay static int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
12765a838052SSatish Balay {
12775a838052SSatish Balay   *bs = 1;
12785a838052SSatish Balay   return 0;
12795a838052SSatish Balay }
12805a838052SSatish Balay 
1281e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
12824dcbc457SBarry Smith {
1283e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
128406763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
12858a047759SSatish Balay   int        start, end, *ai, *aj;
128606763907SSatish Balay   char       *table;
12878a047759SSatish Balay   shift = a->indexshift;
1288e4d965acSSatish Balay   m     = a->m;
1289e4d965acSSatish Balay   ai    = a->i;
12908a047759SSatish Balay   aj    = a->j+shift;
12918a047759SSatish Balay 
1292e4d965acSSatish Balay   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
129306763907SSatish Balay 
129406763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
129506763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
129606763907SSatish Balay 
1297e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1298e4d965acSSatish Balay     /* Initialise the two local arrays */
1299e4d965acSSatish Balay     isz  = 0;
130006763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1301e4d965acSSatish Balay 
1302e4d965acSSatish Balay                 /* Extract the indices, assume there can be duplicate entries */
13034dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
130477c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1305e4d965acSSatish Balay 
1306e4d965acSSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1307e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
130806763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
13094dcbc457SBarry Smith     }
131006763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
131106763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1312e4d965acSSatish Balay 
131304a348a9SBarry Smith     k = 0;
131404a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
131504a348a9SBarry Smith       n = isz;
131606763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1317e4d965acSSatish Balay         row   = nidx[k];
1318e4d965acSSatish Balay         start = ai[row];
1319e4d965acSSatish Balay         end   = ai[row+1];
132004a348a9SBarry Smith         for ( l = start; l<end ; l++){
13218a047759SSatish Balay           val = aj[l] + shift;
132206763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1323e4d965acSSatish Balay         }
1324e4d965acSSatish Balay       }
1325e4d965acSSatish Balay     }
1326537820f0SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1327e4d965acSSatish Balay   }
132804a348a9SBarry Smith   PetscFree(table);
132906763907SSatish Balay   PetscFree(nidx);
1330e4d965acSSatish Balay   return 0;
13314dcbc457SBarry Smith }
133217ab2063SBarry Smith 
1333682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1334682d7d0cSBarry Smith {
1335682d7d0cSBarry Smith   static int called = 0;
1336682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1337682d7d0cSBarry Smith 
1338682d7d0cSBarry Smith   if (called) return 0; else called = 1;
133977c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
134077c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
134177c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
134277c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
134377c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1344682d7d0cSBarry Smith #if defined(HAVE_ESSL)
134577c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1346682d7d0cSBarry Smith #endif
1347682d7d0cSBarry Smith   return 0;
1348682d7d0cSBarry Smith }
1349205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1350682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
135117ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
135217ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1353416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1354416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
135517ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
135617ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
135717ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
135817ab2063SBarry Smith        MatRelax_SeqAIJ,
135917ab2063SBarry Smith        MatTranspose_SeqAIJ,
13607264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1361f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
136217ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
136317ab2063SBarry Smith        MatCompress_SeqAIJ,
136417ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
136517ab2063SBarry Smith        MatGetReordering_SeqAIJ,
136617ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
136717ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
136817ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
136917ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
13703d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1371cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
13727eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1373682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1374f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
13755a838052SSatish Balay        MatScale_SeqAIJ,0,0,
13765a838052SSatish Balay        MatILUDTFactor_SeqAIJ,MatGetBlockSize_SeqAIJ};
137717ab2063SBarry Smith 
137817ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
137917ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
138017ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
138117ab2063SBarry Smith 
138217ab2063SBarry Smith /*@C
1383682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
13840d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
13856e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
13862bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
13872bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
138817ab2063SBarry Smith 
138917ab2063SBarry Smith    Input Parameters:
139017ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
139117ab2063SBarry Smith .  m - number of rows
139217ab2063SBarry Smith .  n - number of columns
139317ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
13942bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
13952bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
139617ab2063SBarry Smith 
139717ab2063SBarry Smith    Output Parameter:
1398416022c9SBarry Smith .  A - the matrix
139917ab2063SBarry Smith 
140017ab2063SBarry Smith    Notes:
140117ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
140217ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
14030002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
140444cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
140517ab2063SBarry Smith 
140617ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1407a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
14083d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
14093d323bbdSBarry Smith    will get TERRIBLE performance, see the users' manual chapter on
14100d15e28bSLois Curfman McInnes    matrices and the file $(PETSC_DIR)/Performance.
141117ab2063SBarry Smith 
1412682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1413682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1414682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
14156c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
14166c7ebb05SLois Curfman McInnes 
14176c7ebb05SLois Curfman McInnes    Options Database Keys:
14186c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
14196e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
14206e62573dSLois Curfman McInnes $        (max limit=5)
14216e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
14226e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
14236e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
142417ab2063SBarry Smith 
142517ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
142617ab2063SBarry Smith @*/
1427416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
142817ab2063SBarry Smith {
1429416022c9SBarry Smith   Mat        B;
1430416022c9SBarry Smith   Mat_SeqAIJ *b;
143169957df2SSatish Balay   int        i, len, ierr, flg;
1432d5d45c9bSBarry Smith 
1433416022c9SBarry Smith   *A                  = 0;
14340452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1435416022c9SBarry Smith   PLogObjectCreate(B);
14360452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
143744cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1438416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1439416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1440416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1441416022c9SBarry Smith   B->factor           = 0;
1442416022c9SBarry Smith   B->lupivotthreshold = 1.0;
14437a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
144469957df2SSatish Balay                           &flg); CHKERRQ(ierr);
14457a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
14467a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
14477a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1448416022c9SBarry Smith   b->row              = 0;
1449416022c9SBarry Smith   b->col              = 0;
1450416022c9SBarry Smith   b->indexshift       = 0;
1451b810aeb4SBarry Smith   b->reallocs         = 0;
145269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
145369957df2SSatish Balay   if (flg) b->indexshift = -1;
145417ab2063SBarry Smith 
145544cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
145644cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
14570452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1458b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
14597b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
14607b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1461416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
146217ab2063SBarry Smith     nz = nz*m;
146317ab2063SBarry Smith   }
146417ab2063SBarry Smith   else {
146517ab2063SBarry Smith     nz = 0;
1466416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
146717ab2063SBarry Smith   }
146817ab2063SBarry Smith 
146917ab2063SBarry Smith   /* allocate the matrix space */
1470416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
14710452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1472416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1473cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1474416022c9SBarry Smith   b->i  = b->j + nz;
1475416022c9SBarry Smith   b->singlemalloc = 1;
147617ab2063SBarry Smith 
1477416022c9SBarry Smith   b->i[0] = -b->indexshift;
147817ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1479416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
148017ab2063SBarry Smith   }
148117ab2063SBarry Smith 
1482416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
14830452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1484416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1485416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
148617ab2063SBarry Smith 
1487416022c9SBarry Smith   b->nz               = 0;
1488416022c9SBarry Smith   b->maxnz            = nz;
1489416022c9SBarry Smith   b->sorted           = 0;
1490416022c9SBarry Smith   b->roworiented      = 1;
1491416022c9SBarry Smith   b->nonew            = 0;
1492416022c9SBarry Smith   b->diag             = 0;
1493416022c9SBarry Smith   b->solve_work       = 0;
149471bd300dSLois Curfman McInnes   b->spptr            = 0;
1495754ec7b1SSatish Balay   b->inode.node_count = 0;
1496754ec7b1SSatish Balay   b->inode.size       = 0;
14976c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
14986c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
1499*4e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
150017ab2063SBarry Smith 
1501416022c9SBarry Smith   *A = B;
1502*4e220ebcSLois Curfman McInnes 
15034b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
15044b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
150569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
150669957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
15074b14c69eSBarry Smith #endif
150869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
150969957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
151069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
151169957df2SSatish Balay   if (flg) {
1512416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1513416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
151417ab2063SBarry Smith   }
151569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
151669957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
151717ab2063SBarry Smith   return 0;
151817ab2063SBarry Smith }
151917ab2063SBarry Smith 
15203d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
152117ab2063SBarry Smith {
1522416022c9SBarry Smith   Mat        C;
1523416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
152408480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
152517ab2063SBarry Smith 
15264043dd9cSLois Curfman McInnes   *B = 0;
15270452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1528416022c9SBarry Smith   PLogObjectCreate(C);
15290452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
153041c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1531416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1532416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1533416022c9SBarry Smith   C->factor     = A->factor;
1534416022c9SBarry Smith   c->row        = 0;
1535416022c9SBarry Smith   c->col        = 0;
1536416022c9SBarry Smith   c->indexshift = shift;
1537c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
153817ab2063SBarry Smith 
153944cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
154044cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
154144cd7ae7SLois Curfman McInnes   C->M          = a->m;
154244cd7ae7SLois Curfman McInnes   C->N          = a->n;
154317ab2063SBarry Smith 
15440452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
15450452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
154617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1547416022c9SBarry Smith     c->imax[i] = a->imax[i];
1548416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
154917ab2063SBarry Smith   }
155017ab2063SBarry Smith 
155117ab2063SBarry Smith   /* allocate the matrix space */
1552416022c9SBarry Smith   c->singlemalloc = 1;
1553416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
15540452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1555416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1556416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1557416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
155817ab2063SBarry Smith   if (m > 0) {
1559416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
156008480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1561416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
156217ab2063SBarry Smith     }
156308480c60SBarry Smith   }
156417ab2063SBarry Smith 
1565416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1566416022c9SBarry Smith   c->sorted      = a->sorted;
1567416022c9SBarry Smith   c->roworiented = a->roworiented;
1568416022c9SBarry Smith   c->nonew       = a->nonew;
15697a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
157017ab2063SBarry Smith 
1571416022c9SBarry Smith   if (a->diag) {
15720452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1573416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
157417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1575416022c9SBarry Smith       c->diag[i] = a->diag[i];
157617ab2063SBarry Smith     }
157717ab2063SBarry Smith   }
1578416022c9SBarry Smith   else c->diag          = 0;
15796c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
15806c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1581754ec7b1SSatish Balay   if (a->inode.size){
1582754ec7b1SSatish Balay     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1583754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1584754ec7b1SSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1585754ec7b1SSatish Balay   } else {
1586754ec7b1SSatish Balay     c->inode.size       = 0;
1587754ec7b1SSatish Balay     c->inode.node_count = 0;
1588754ec7b1SSatish Balay   }
1589416022c9SBarry Smith   c->nz                 = a->nz;
1590416022c9SBarry Smith   c->maxnz              = a->maxnz;
1591416022c9SBarry Smith   c->solve_work         = 0;
159276dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1593754ec7b1SSatish Balay 
1594416022c9SBarry Smith   *B = C;
159517ab2063SBarry Smith   return 0;
159617ab2063SBarry Smith }
159717ab2063SBarry Smith 
159819bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
159917ab2063SBarry Smith {
1600416022c9SBarry Smith   Mat_SeqAIJ   *a;
1601416022c9SBarry Smith   Mat          B;
160217699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1603bcd2baecSBarry Smith   MPI_Comm     comm;
160417ab2063SBarry Smith 
160519bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
160617699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
160717699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
160890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
160977c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
161048d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
161117ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
161217ab2063SBarry Smith 
161317ab2063SBarry Smith   /* read in row lengths */
16140452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
161577c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
161617ab2063SBarry Smith 
161717ab2063SBarry Smith   /* create our matrix */
1618416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1619416022c9SBarry Smith   B = *A;
1620416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1621416022c9SBarry Smith   shift = a->indexshift;
162217ab2063SBarry Smith 
162317ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
162477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
162517ab2063SBarry Smith   if (shift) {
162617ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1627416022c9SBarry Smith       a->j[i] += 1;
162817ab2063SBarry Smith     }
162917ab2063SBarry Smith   }
163017ab2063SBarry Smith 
163117ab2063SBarry Smith   /* read in nonzero values */
163277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
163317ab2063SBarry Smith 
163417ab2063SBarry Smith   /* set matrix "i" values */
1635416022c9SBarry Smith   a->i[0] = -shift;
163617ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1637416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1638416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
163917ab2063SBarry Smith   }
16400452661fSBarry Smith   PetscFree(rowlengths);
164117ab2063SBarry Smith 
16426d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
16436d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
164417ab2063SBarry Smith   return 0;
164517ab2063SBarry Smith }
164617ab2063SBarry Smith 
1647686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
16487264ac53SSatish Balay {
16497264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
16507264ac53SSatish Balay 
1651bcd2baecSBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type");
16527264ac53SSatish Balay 
16537264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
16547264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1655bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
165677c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1657bcd2baecSBarry Smith   }
16587264ac53SSatish Balay 
16597264ac53SSatish Balay   /* if the a->i are the same */
1660bcd2baecSBarry Smith   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
166177c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
16627264ac53SSatish Balay   }
16637264ac53SSatish Balay 
16647264ac53SSatish Balay   /* if a->j are the same */
1665bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
166677c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1667bcd2baecSBarry Smith   }
1668bcd2baecSBarry Smith 
1669bcd2baecSBarry Smith   /* if a->a are the same */
167019bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
167177c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
16727264ac53SSatish Balay   }
167377c4ece6SBarry Smith   *flg = PETSC_TRUE;
16747264ac53SSatish Balay   return 0;
16757264ac53SSatish Balay 
16767264ac53SSatish Balay }
1677