xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 9b1297e1131a6dc2f27c05e188d7eebf48061259)
16d84be18SBarry Smith 
26945ee14SBarry Smith 
317ab2063SBarry Smith #ifndef lint
4*9b1297e1SSatish Balay static char vcid[] = "$Id: aij.c,v 1.192 1996/11/01 16:47:14 balay Exp balay $";
517ab2063SBarry Smith #endif
617ab2063SBarry Smith 
7d5d45c9bSBarry Smith /*
85a838052SSatish Balay B    Defines the basic matrix operations for the AIJ (compressed row)
9d5d45c9bSBarry Smith   matrix storage format.
10d5d45c9bSBarry Smith */
1170f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
12f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
13f5eb4b81SSatish Balay #include "src/inline/spops.h"
14e4d965acSSatish Balay #include "petsc.h"
15f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1617ab2063SBarry Smith 
17a2ce50c7SBarry Smith /*
18a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
19a2ce50c7SBarry Smith */
20a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
21a2ce50c7SBarry Smith {
22a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
23a2ce50c7SBarry Smith   int        ierr = 1;
24a2ce50c7SBarry Smith 
25a2ce50c7SBarry Smith   SETERRQ(ierr,"MatILUDTFactor_SeqAIJ:Not implemented");
26a2ce50c7SBarry Smith }
27a2ce50c7SBarry Smith 
28bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
2917ab2063SBarry Smith 
3031625ec5SSatish Balay static int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
316945ee14SBarry Smith                            PetscTruth *done)
3217ab2063SBarry Smith {
33416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
346945ee14SBarry Smith   int        ierr,i,ishift;
3517ab2063SBarry Smith 
3631625ec5SSatish Balay   *m     = A->m;
376945ee14SBarry Smith   if (!ia) return 0;
386945ee14SBarry Smith   ishift = a->indexshift;
396945ee14SBarry Smith   if (symmetric) {
4031625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
416945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
4231625ec5SSatish Balay     int nz = a->i[a->m];
433b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
4431625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
453b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
463b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
4731625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
486945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
4931625ec5SSatish Balay     int nz = a->i[a->m] + 1;
503b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
5131625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
523b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
533b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
5431625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
556945ee14SBarry Smith   } else {
566945ee14SBarry Smith     *ia = a->i; *ja = a->j;
57a2ce50c7SBarry Smith   }
58a2ce50c7SBarry Smith 
59a2744918SBarry Smith   return 0;
60a2744918SBarry Smith }
61a2744918SBarry Smith 
623b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
636945ee14SBarry Smith                                PetscTruth *done)
646945ee14SBarry Smith {
656945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
663b2fbd54SBarry Smith   int        ishift = a->indexshift;
676945ee14SBarry Smith 
686945ee14SBarry Smith   if (!ia) return 0;
693b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
706945ee14SBarry Smith     PetscFree(*ia);
716945ee14SBarry Smith     PetscFree(*ja);
72bcd2baecSBarry Smith   }
7317ab2063SBarry Smith   return 0;
7417ab2063SBarry Smith }
7517ab2063SBarry Smith 
763b2fbd54SBarry Smith static int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
773b2fbd54SBarry Smith                            PetscTruth *done)
783b2fbd54SBarry Smith {
793b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
80a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
81a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
823b2fbd54SBarry Smith 
833b2fbd54SBarry Smith   *nn     = A->n;
843b2fbd54SBarry Smith   if (!ia) return 0;
853b2fbd54SBarry Smith   if (symmetric) {
86179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
873b2fbd54SBarry Smith   } else {
8861d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
893b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
903b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
91a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
923b2fbd54SBarry Smith     jj = a->j;
933b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
943b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
953b2fbd54SBarry Smith     }
963b2fbd54SBarry Smith     cia[0] = oshift;
973b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
983b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
993b2fbd54SBarry Smith     }
1003b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1013b2fbd54SBarry Smith     jj = a->j;
102a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
103a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
104a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1053b2fbd54SBarry Smith         col = *jj++ + ishift;
1063b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1073b2fbd54SBarry Smith       }
1083b2fbd54SBarry Smith     }
1093b2fbd54SBarry Smith     PetscFree(collengths);
1103b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1113b2fbd54SBarry Smith   }
1123b2fbd54SBarry Smith 
1133b2fbd54SBarry Smith   return 0;
1143b2fbd54SBarry Smith }
1153b2fbd54SBarry Smith 
1163b2fbd54SBarry Smith static int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1173b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1183b2fbd54SBarry Smith {
1193b2fbd54SBarry Smith   if (!ia) return 0;
1203b2fbd54SBarry Smith 
1213b2fbd54SBarry Smith   PetscFree(*ia);
1223b2fbd54SBarry Smith   PetscFree(*ja);
1233b2fbd54SBarry Smith 
1243b2fbd54SBarry Smith   return 0;
1253b2fbd54SBarry Smith }
1263b2fbd54SBarry Smith 
127227d817aSBarry Smith #define CHUNKSIZE   15
12817ab2063SBarry Smith 
12961d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
13017ab2063SBarry Smith {
131416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
132416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1334b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
134d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
135416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
13617ab2063SBarry Smith 
13717ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
138416022c9SBarry Smith     row  = im[k];
1393b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
14017ab2063SBarry Smith     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
141416022c9SBarry Smith     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
1423b2fbd54SBarry Smith #endif
14317ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
14417ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
145416022c9SBarry Smith     low = 0;
14617ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1473b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
148416022c9SBarry Smith       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
149416022c9SBarry Smith       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
1503b2fbd54SBarry Smith #endif
1514b0e389bSBarry Smith       col = in[l] - shift;
1524b0e389bSBarry Smith       if (roworiented) {
1534b0e389bSBarry Smith         value = *v++;
1544b0e389bSBarry Smith       }
1554b0e389bSBarry Smith       else {
1564b0e389bSBarry Smith         value = v[k + l*m];
1574b0e389bSBarry Smith       }
158416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
159416022c9SBarry Smith       while (high-low > 5) {
160416022c9SBarry Smith         t = (low+high)/2;
161416022c9SBarry Smith         if (rp[t] > col) high = t;
162416022c9SBarry Smith         else             low  = t;
16317ab2063SBarry Smith       }
164416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
16517ab2063SBarry Smith         if (rp[i] > col) break;
16617ab2063SBarry Smith         if (rp[i] == col) {
167416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
16817ab2063SBarry Smith           else                  ap[i] = value;
16917ab2063SBarry Smith           goto noinsert;
17017ab2063SBarry Smith         }
17117ab2063SBarry Smith       }
17217ab2063SBarry Smith       if (nonew) goto noinsert;
17317ab2063SBarry Smith       if (nrow >= rmax) {
17417ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
175416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
17617ab2063SBarry Smith         Scalar *new_a;
17717ab2063SBarry Smith 
17817ab2063SBarry Smith         /* malloc new storage space */
179416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
1800452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
18117ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
18217ab2063SBarry Smith         new_i   = new_j + new_nz;
18317ab2063SBarry Smith 
18417ab2063SBarry Smith         /* copy over old data into new slots */
18517ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
186416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
187416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
188416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
189416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
19017ab2063SBarry Smith                                                            len*sizeof(int));
191416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
192416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
19317ab2063SBarry Smith                                                            len*sizeof(Scalar));
19417ab2063SBarry Smith         /* free up old matrix storage */
1950452661fSBarry Smith         PetscFree(a->a);
1960452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
197416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
198416022c9SBarry Smith         a->singlemalloc = 1;
19917ab2063SBarry Smith 
20017ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
201416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
202416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
203416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
204b810aeb4SBarry Smith         a->reallocs++;
20517ab2063SBarry Smith       }
206416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
207416022c9SBarry Smith       /* shift up all the later entries in this row */
208416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
20917ab2063SBarry Smith         rp[ii+1] = rp[ii];
21017ab2063SBarry Smith         ap[ii+1] = ap[ii];
21117ab2063SBarry Smith       }
21217ab2063SBarry Smith       rp[i] = col;
21317ab2063SBarry Smith       ap[i] = value;
21417ab2063SBarry Smith       noinsert:;
215416022c9SBarry Smith       low = i + 1;
21617ab2063SBarry Smith     }
21717ab2063SBarry Smith     ailen[row] = nrow;
21817ab2063SBarry Smith   }
21917ab2063SBarry Smith   return 0;
22017ab2063SBarry Smith }
22117ab2063SBarry Smith 
2227eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2237eb43aa7SLois Curfman McInnes {
2247eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
225b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2267eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2277eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2287eb43aa7SLois Curfman McInnes 
2297eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2307eb43aa7SLois Curfman McInnes     row  = im[k];
2317eb43aa7SLois Curfman McInnes     if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row");
2327eb43aa7SLois Curfman McInnes     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large");
2337eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2347eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2357eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
2367eb43aa7SLois Curfman McInnes       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column");
2377eb43aa7SLois Curfman McInnes       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large");
2387eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2397eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2407eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2417eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2427eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2437eb43aa7SLois Curfman McInnes         else             low  = t;
2447eb43aa7SLois Curfman McInnes       }
2457eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2467eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2477eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
248b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2497eb43aa7SLois Curfman McInnes           goto finished;
2507eb43aa7SLois Curfman McInnes         }
2517eb43aa7SLois Curfman McInnes       }
252b49de8d1SLois Curfman McInnes       *v++ = zero;
2537eb43aa7SLois Curfman McInnes       finished:;
2547eb43aa7SLois Curfman McInnes     }
2557eb43aa7SLois Curfman McInnes   }
2567eb43aa7SLois Curfman McInnes   return 0;
2577eb43aa7SLois Curfman McInnes }
2587eb43aa7SLois Curfman McInnes 
25917ab2063SBarry Smith #include "draw.h"
26017ab2063SBarry Smith #include "pinclude/pviewer.h"
26177c4ece6SBarry Smith #include "sys.h"
26217ab2063SBarry Smith 
263416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
26417ab2063SBarry Smith {
265416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
266416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
26717ab2063SBarry Smith 
26890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2690452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
270416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
271416022c9SBarry Smith   col_lens[1] = a->m;
272416022c9SBarry Smith   col_lens[2] = a->n;
273416022c9SBarry Smith   col_lens[3] = a->nz;
274416022c9SBarry Smith 
275416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
276416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
277416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
27817ab2063SBarry Smith   }
27977c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
2800452661fSBarry Smith   PetscFree(col_lens);
281416022c9SBarry Smith 
282416022c9SBarry Smith   /* store column indices (zero start index) */
283416022c9SBarry Smith   if (a->indexshift) {
284416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
28517ab2063SBarry Smith   }
28677c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr);
287416022c9SBarry Smith   if (a->indexshift) {
288416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
28917ab2063SBarry Smith   }
290416022c9SBarry Smith 
291416022c9SBarry Smith   /* store nonzero values */
29277c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
29317ab2063SBarry Smith   return 0;
29417ab2063SBarry Smith }
295416022c9SBarry Smith 
296416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
297416022c9SBarry Smith {
298416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
299496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
30017ab2063SBarry Smith   FILE        *fd;
30117ab2063SBarry Smith   char        *outputname;
30217ab2063SBarry Smith 
30390ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
304416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
30590ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
306a93ec695SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
30795e01e2fSLois Curfman McInnes     return 0;
30895e01e2fSLois Curfman McInnes   }
309a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
310496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
311496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
312496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
31395e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
31495e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
31517ab2063SBarry Smith   }
316a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
317416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3184e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
3194e220ebcSLois Curfman McInnes     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz);
32017ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
32117ab2063SBarry Smith 
32217ab2063SBarry Smith     for (i=0; i<m; i++) {
323416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
32417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
3256945ee14SBarry Smith         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]),
326416022c9SBarry Smith                    imag(a->a[j]));
32717ab2063SBarry Smith #else
3287a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
32917ab2063SBarry Smith #endif
33017ab2063SBarry Smith       }
33117ab2063SBarry Smith     }
33217ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
33317ab2063SBarry Smith   }
334a93ec695SBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
33544cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
33644cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
33744cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
33844cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
33944cd7ae7SLois Curfman McInnes         if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0)
34044cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
34144cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
34244cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
34344cd7ae7SLois Curfman McInnes #else
34444cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
34544cd7ae7SLois Curfman McInnes #endif
34644cd7ae7SLois Curfman McInnes       }
34744cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
34844cd7ae7SLois Curfman McInnes     }
34944cd7ae7SLois Curfman McInnes   }
35017ab2063SBarry Smith   else {
35117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
35217ab2063SBarry Smith       fprintf(fd,"row %d:",i);
353416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
35417ab2063SBarry Smith #if defined(PETSC_COMPLEX)
355416022c9SBarry Smith         if (imag(a->a[j]) != 0.0) {
356416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
35717ab2063SBarry Smith         }
35817ab2063SBarry Smith         else {
359416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
36017ab2063SBarry Smith         }
36117ab2063SBarry Smith #else
362416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
36317ab2063SBarry Smith #endif
36417ab2063SBarry Smith       }
36517ab2063SBarry Smith       fprintf(fd,"\n");
36617ab2063SBarry Smith     }
36717ab2063SBarry Smith   }
36817ab2063SBarry Smith   fflush(fd);
369416022c9SBarry Smith   return 0;
370416022c9SBarry Smith }
371416022c9SBarry Smith 
372416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
373416022c9SBarry Smith {
374416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
375cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
376cddf8d76SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
377bcd2baecSBarry Smith   Draw        draw;
378cddf8d76SBarry Smith   DrawButton  button;
37919bcc07fSBarry Smith   PetscTruth  isnull;
380cddf8d76SBarry Smith 
381bcd2baecSBarry Smith   ViewerDrawGetDraw(viewer,&draw);
38219bcc07fSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
38319bcc07fSBarry Smith 
384416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
385416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
386416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
387416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
388cddf8d76SBarry Smith   color = DRAW_BLUE;
389416022c9SBarry Smith   for ( i=0; i<m; i++ ) {
390cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
391416022c9SBarry 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);
399cddf8d76SBarry Smith     }
400cddf8d76SBarry Smith   }
401cddf8d76SBarry Smith   color = DRAW_CYAN;
402cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
403cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
404cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
405cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
406cddf8d76SBarry Smith       if (a->a[j] !=  0.) continue;
407cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
408cddf8d76SBarry Smith     }
409cddf8d76SBarry Smith   }
410cddf8d76SBarry Smith   color = DRAW_RED;
411cddf8d76SBarry Smith   for ( i=0; i<m; i++ ) {
412cddf8d76SBarry Smith     y_l = m - i - 1.0; y_r = y_l + 1.0;
413cddf8d76SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
414cddf8d76SBarry Smith       x_l = a->j[j] + shift; x_r = x_l + 1.0;
415cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
416cddf8d76SBarry Smith       if (real(a->a[j]) <=  0.) continue;
417cddf8d76SBarry Smith #else
418cddf8d76SBarry Smith       if (a->a[j] <=  0.) continue;
419cddf8d76SBarry Smith #endif
420cddf8d76SBarry Smith       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
421416022c9SBarry Smith     }
422416022c9SBarry Smith   }
423416022c9SBarry Smith   DrawFlush(draw);
424cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
425cddf8d76SBarry Smith   if (pause >= 0) { PetscSleep(pause); return 0;}
426cddf8d76SBarry Smith 
427cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
4286945ee14SBarry Smith   ierr = DrawCheckResizedWindow(draw);
429cddf8d76SBarry Smith   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
430cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
431cddf8d76SBarry Smith     DrawClear(draw);
432cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
433cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
434cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
435cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
436cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
437cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
438cddf8d76SBarry Smith     w *= scale; h *= scale;
439cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
440cddf8d76SBarry Smith     color = DRAW_BLUE;
441cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
442cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
443cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
444cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
445cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
446cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
447cddf8d76SBarry Smith #else
448cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
449cddf8d76SBarry Smith #endif
450cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
451cddf8d76SBarry Smith       }
452cddf8d76SBarry Smith     }
453cddf8d76SBarry Smith     color = DRAW_CYAN;
454cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
455cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
456cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
457cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
458cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
459cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
460cddf8d76SBarry Smith       }
461cddf8d76SBarry Smith     }
462cddf8d76SBarry Smith     color = DRAW_RED;
463cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
464cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
465cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
466cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
467cddf8d76SBarry Smith #if defined(PETSC_COMPLEX)
468cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
469cddf8d76SBarry Smith #else
470cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
471cddf8d76SBarry Smith #endif
472cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
473cddf8d76SBarry Smith       }
474cddf8d76SBarry Smith     }
4756945ee14SBarry Smith     ierr = DrawCheckResizedWindow(draw);
476cddf8d76SBarry Smith     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
477cddf8d76SBarry Smith   }
478416022c9SBarry Smith   return 0;
479416022c9SBarry Smith }
480416022c9SBarry Smith 
481416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
482416022c9SBarry Smith {
483416022c9SBarry Smith   Mat         A = (Mat) obj;
484416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
485bcd2baecSBarry Smith   ViewerType  vtype;
486bcd2baecSBarry Smith   int         ierr;
487416022c9SBarry Smith 
488bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
489bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
490416022c9SBarry Smith     return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
491416022c9SBarry Smith   }
492bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
493416022c9SBarry Smith     return MatView_SeqAIJ_ASCII(A,viewer);
494416022c9SBarry Smith   }
495bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
496416022c9SBarry Smith     return MatView_SeqAIJ_Binary(A,viewer);
497416022c9SBarry Smith   }
498bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
499bcd2baecSBarry Smith     return MatView_SeqAIJ_Draw(A,viewer);
50017ab2063SBarry Smith   }
50117ab2063SBarry Smith   return 0;
50217ab2063SBarry Smith }
50319bcc07fSBarry Smith 
504c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
505416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
50617ab2063SBarry Smith {
507416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
50841c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
509416022c9SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
510416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
51117ab2063SBarry Smith 
5126d4a8577SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
51317ab2063SBarry Smith 
51417ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
515416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
51617ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
51717ab2063SBarry Smith     if (fshift) {
518416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
51917ab2063SBarry Smith       N = ailen[i];
52017ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
52117ab2063SBarry Smith         ip[j-fshift] = ip[j];
52217ab2063SBarry Smith         ap[j-fshift] = ap[j];
52317ab2063SBarry Smith       }
52417ab2063SBarry Smith     }
52517ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
52617ab2063SBarry Smith   }
52717ab2063SBarry Smith   if (m) {
52817ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
52917ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
53017ab2063SBarry Smith   }
53117ab2063SBarry Smith   /* reset ilen and imax for each row */
53217ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
53317ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
53417ab2063SBarry Smith   }
535416022c9SBarry Smith   a->nz = ai[m] + shift;
53617ab2063SBarry Smith 
53717ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
538416022c9SBarry Smith   if (fshift && a->diag) {
5390452661fSBarry Smith     PetscFree(a->diag);
540416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
541416022c9SBarry Smith     a->diag = 0;
54217ab2063SBarry Smith   }
5434e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
5444e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
5454e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
546b810aeb4SBarry Smith            a->reallocs);
5474e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
5484e220ebcSLois Curfman McInnes 
54976dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
55041c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
55117ab2063SBarry Smith   return 0;
55217ab2063SBarry Smith }
55317ab2063SBarry Smith 
554416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A)
55517ab2063SBarry Smith {
556416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
557cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
55817ab2063SBarry Smith   return 0;
55917ab2063SBarry Smith }
560416022c9SBarry Smith 
56117ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj)
56217ab2063SBarry Smith {
563416022c9SBarry Smith   Mat        A  = (Mat) obj;
564416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
565d5d45c9bSBarry Smith 
56617ab2063SBarry Smith #if defined(PETSC_LOG)
567416022c9SBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
56817ab2063SBarry Smith #endif
5690452661fSBarry Smith   PetscFree(a->a);
5700452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
5710452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
5720452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
5730452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
5740452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
57576dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
5760452661fSBarry Smith   PetscFree(a);
577f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
578f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
57917ab2063SBarry Smith   return 0;
58017ab2063SBarry Smith }
58117ab2063SBarry Smith 
582416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A)
58317ab2063SBarry Smith {
58417ab2063SBarry Smith   return 0;
58517ab2063SBarry Smith }
58617ab2063SBarry Smith 
587416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op)
58817ab2063SBarry Smith {
589416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
5906d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)              a->roworiented = 1;
5916d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)           a->roworiented = 0;
5926d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)            a->sorted      = 1;
5936d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
5946d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
5956d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
5966d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
5976d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
5986d4a8577SBarry Smith            op == MAT_YES_NEW_DIAGONALS)
59994a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
6006d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
6016d4a8577SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:MAT_NO_NEW_DIAGONALS");}
6026d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
6036d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
6046d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
6056d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
6066d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
607e2f28af5SBarry Smith   else
6084d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
60917ab2063SBarry Smith   return 0;
61017ab2063SBarry Smith }
61117ab2063SBarry Smith 
612416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
61317ab2063SBarry Smith {
614416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
615416022c9SBarry Smith   int        i,j, n,shift = a->indexshift;
61617ab2063SBarry Smith   Scalar     *x, zero = 0.0;
61717ab2063SBarry Smith 
61817ab2063SBarry Smith   VecSet(&zero,v);
61917ab2063SBarry Smith   VecGetArray(v,&x); VecGetLocalSize(v,&n);
620416022c9SBarry Smith   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
621416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
622416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
623416022c9SBarry Smith       if (a->j[j]+shift == i) {
624416022c9SBarry Smith         x[i] = a->a[j];
62517ab2063SBarry Smith         break;
62617ab2063SBarry Smith       }
62717ab2063SBarry Smith     }
62817ab2063SBarry Smith   }
62917ab2063SBarry Smith   return 0;
63017ab2063SBarry Smith }
63117ab2063SBarry Smith 
63217ab2063SBarry Smith /* -------------------------------------------------------*/
63317ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
63417ab2063SBarry Smith /* -------------------------------------------------------*/
63544cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
63617ab2063SBarry Smith {
637416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
63817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
639416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift;
64017ab2063SBarry Smith 
64117ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
642cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
64317ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
64417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
645416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
646416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
647416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
64817ab2063SBarry Smith     alpha = x[i];
64917ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
65017ab2063SBarry Smith   }
651416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
65217ab2063SBarry Smith   return 0;
65317ab2063SBarry Smith }
654d5d45c9bSBarry Smith 
65544cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
65617ab2063SBarry Smith {
657416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
65817ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
659416022c9SBarry Smith   int        m = a->m, n, i, *idx,shift = a->indexshift;
66017ab2063SBarry Smith 
66117ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
66217ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
66317ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
66417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
665416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
666416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
667416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
66817ab2063SBarry Smith     alpha = x[i];
66917ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
67017ab2063SBarry Smith   }
67117ab2063SBarry Smith   return 0;
67217ab2063SBarry Smith }
67317ab2063SBarry Smith 
67444cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
67517ab2063SBarry Smith {
676416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
67717ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
678416022c9SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
67917ab2063SBarry Smith 
68017ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
68117ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
682416022c9SBarry Smith   idx  = a->j;
683416022c9SBarry Smith   v    = a->a;
684416022c9SBarry Smith   ii   = a->i;
68517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
686416022c9SBarry Smith     n    = ii[1] - ii[0]; ii++;
68717ab2063SBarry Smith     sum  = 0.0;
68817ab2063SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
68917ab2063SBarry Smith     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
690416022c9SBarry Smith     while (n--) sum += *v++ * x[*idx++];
69117ab2063SBarry Smith     y[i] = sum;
69217ab2063SBarry Smith   }
693416022c9SBarry Smith   PLogFlops(2*a->nz - m);
69417ab2063SBarry Smith   return 0;
69517ab2063SBarry Smith }
69617ab2063SBarry Smith 
69744cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
69817ab2063SBarry Smith {
699416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
70017ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
701cddf8d76SBarry Smith   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
70217ab2063SBarry Smith 
70317ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
70417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
705cddf8d76SBarry Smith   idx  = a->j;
706cddf8d76SBarry Smith   v    = a->a;
707cddf8d76SBarry Smith   ii   = a->i;
70817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
709cddf8d76SBarry Smith     n    = ii[1] - ii[0]; ii++;
71017ab2063SBarry Smith     sum  = y[i];
711cddf8d76SBarry Smith     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
712cddf8d76SBarry Smith     while (n--) sum += *v++ * x[*idx++];
71317ab2063SBarry Smith     z[i] = sum;
71417ab2063SBarry Smith   }
715416022c9SBarry Smith   PLogFlops(2*a->nz);
71617ab2063SBarry Smith   return 0;
71717ab2063SBarry Smith }
71817ab2063SBarry Smith 
71917ab2063SBarry Smith /*
72017ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
72117ab2063SBarry Smith */
72217ab2063SBarry Smith 
723416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
72417ab2063SBarry Smith {
725416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
726416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
72717ab2063SBarry Smith 
7280452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
729416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
730416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
731416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
732416022c9SBarry Smith       if (a->j[j]+shift == i) {
73317ab2063SBarry Smith         diag[i] = j - shift;
73417ab2063SBarry Smith         break;
73517ab2063SBarry Smith       }
73617ab2063SBarry Smith     }
73717ab2063SBarry Smith   }
738416022c9SBarry Smith   a->diag = diag;
73917ab2063SBarry Smith   return 0;
74017ab2063SBarry Smith }
74117ab2063SBarry Smith 
74244cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
74317ab2063SBarry Smith                            double fshift,int its,Vec xx)
74417ab2063SBarry Smith {
745416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
746416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
747d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
74817ab2063SBarry Smith 
74917ab2063SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
750416022c9SBarry Smith   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
751416022c9SBarry Smith   diag = a->diag;
752416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
75317ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
75417ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
75517ab2063SBarry Smith     bs = b + shift;
75617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
757416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
758416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
759416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
760416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
76117ab2063SBarry Smith         sum  = b[i]*d/omega;
76217ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
76317ab2063SBarry Smith         x[i] = sum;
76417ab2063SBarry Smith     }
76517ab2063SBarry Smith     return 0;
76617ab2063SBarry Smith   }
76717ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
76817ab2063SBarry Smith     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
76917ab2063SBarry Smith   }
770416022c9SBarry Smith   else if (flag & SOR_EISENSTAT) {
77117ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
77217ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
77317ab2063SBarry Smith 
77417ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
77517ab2063SBarry Smith 
77617ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
77717ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
77817ab2063SBarry Smith     is the relaxation factor.
77917ab2063SBarry Smith     */
7800452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
78117ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
78217ab2063SBarry Smith 
78317ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
78417ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
785416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
786416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
787416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
788416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
78917ab2063SBarry Smith       sum  = b[i];
79017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
79117ab2063SBarry Smith       x[i] = omega*(sum/d);
79217ab2063SBarry Smith     }
79317ab2063SBarry Smith 
79417ab2063SBarry Smith     /*  t = b - (2*E - D)x */
795416022c9SBarry Smith     v = a->a;
79617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
79717ab2063SBarry Smith 
79817ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
799416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
800416022c9SBarry Smith     diag = a->diag;
80117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
802416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
803416022c9SBarry Smith       n    = diag[i] - a->i[i];
804416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
805416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
80617ab2063SBarry Smith       sum  = t[i];
80717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
80817ab2063SBarry Smith       t[i] = omega*(sum/d);
80917ab2063SBarry Smith     }
81017ab2063SBarry Smith 
81117ab2063SBarry Smith     /*  x = x + t */
81217ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
8130452661fSBarry Smith     PetscFree(t);
81417ab2063SBarry Smith     return 0;
81517ab2063SBarry Smith   }
81617ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
81717ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
81817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
819416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
820416022c9SBarry Smith         n    = diag[i] - a->i[i];
821416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
822416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
82317ab2063SBarry Smith         sum  = b[i];
82417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
82517ab2063SBarry Smith         x[i] = omega*(sum/d);
82617ab2063SBarry Smith       }
82717ab2063SBarry Smith       xb = x;
82817ab2063SBarry Smith     }
82917ab2063SBarry Smith     else xb = b;
83017ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
83117ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
83217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
833416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
83417ab2063SBarry Smith       }
83517ab2063SBarry Smith     }
83617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
83717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
838416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
839416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
840416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
841416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
84217ab2063SBarry Smith         sum  = xb[i];
84317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
84417ab2063SBarry Smith         x[i] = omega*(sum/d);
84517ab2063SBarry Smith       }
84617ab2063SBarry Smith     }
84717ab2063SBarry Smith     its--;
84817ab2063SBarry Smith   }
84917ab2063SBarry Smith   while (its--) {
85017ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
85117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
852416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
853416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
854416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
855416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
85617ab2063SBarry Smith         sum  = b[i];
85717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
8587e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
85917ab2063SBarry Smith       }
86017ab2063SBarry Smith     }
86117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
86217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
863416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
864416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
865416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
866416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
86717ab2063SBarry Smith         sum  = b[i];
86817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
8697e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
87017ab2063SBarry Smith       }
87117ab2063SBarry Smith     }
87217ab2063SBarry Smith   }
87317ab2063SBarry Smith   return 0;
87417ab2063SBarry Smith }
87517ab2063SBarry Smith 
8764e220ebcSLois Curfman McInnes static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
87717ab2063SBarry Smith {
878416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
8794e220ebcSLois Curfman McInnes 
8804e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
8814e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
8824e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
8834e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
8844e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
8854e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
8864e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
8874e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
8884e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
8894e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
8904e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
8914e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
8924e220ebcSLois Curfman McInnes   info->memory         = A->mem;
8934e220ebcSLois Curfman McInnes   if (A->factor) {
8944e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
8954e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
8964e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
8974e220ebcSLois Curfman McInnes   } else {
8984e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
8994e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
9004e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
9014e220ebcSLois Curfman McInnes   }
90217ab2063SBarry Smith   return 0;
90317ab2063SBarry Smith }
90417ab2063SBarry Smith 
90517ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
90617ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
90717ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
90817ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
90917ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
91017ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
91117ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
91217ab2063SBarry Smith 
91317ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
91417ab2063SBarry Smith {
915416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
916416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
91717ab2063SBarry Smith 
91877c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
91917ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
92017ab2063SBarry Smith   if (diag) {
92117ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
922416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
923416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
924416022c9SBarry Smith         a->ilen[rows[i]] = 1;
925416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
926416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
92717ab2063SBarry Smith       }
92817ab2063SBarry Smith       else {
92917ab2063SBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
93017ab2063SBarry Smith         CHKERRQ(ierr);
93117ab2063SBarry Smith       }
93217ab2063SBarry Smith     }
93317ab2063SBarry Smith   }
93417ab2063SBarry Smith   else {
93517ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
936416022c9SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
937416022c9SBarry Smith       a->ilen[rows[i]] = 0;
93817ab2063SBarry Smith     }
93917ab2063SBarry Smith   }
94017ab2063SBarry Smith   ISRestoreIndices(is,&rows);
9416d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9426d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
94317ab2063SBarry Smith   return 0;
94417ab2063SBarry Smith }
94517ab2063SBarry Smith 
946416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
94717ab2063SBarry Smith {
948416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
949416022c9SBarry Smith   *m = a->m; *n = a->n;
95017ab2063SBarry Smith   return 0;
95117ab2063SBarry Smith }
95217ab2063SBarry Smith 
953416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
95417ab2063SBarry Smith {
955416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
956416022c9SBarry Smith   *m = 0; *n = a->m;
95717ab2063SBarry Smith   return 0;
95817ab2063SBarry Smith }
9594e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
96017ab2063SBarry Smith {
961416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
962c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
96317ab2063SBarry Smith 
964416022c9SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
96517ab2063SBarry Smith 
966416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
967416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
96817ab2063SBarry Smith   if (idx) {
969416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
9704e093b46SBarry Smith     if (*nz && shift) {
9710452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
97217ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
9734e093b46SBarry Smith     } else if (*nz) {
9744e093b46SBarry Smith       *idx = itmp;
97517ab2063SBarry Smith     }
97617ab2063SBarry Smith     else *idx = 0;
97717ab2063SBarry Smith   }
97817ab2063SBarry Smith   return 0;
97917ab2063SBarry Smith }
98017ab2063SBarry Smith 
9814e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
98217ab2063SBarry Smith {
9834e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
9844e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
98517ab2063SBarry Smith   return 0;
98617ab2063SBarry Smith }
98717ab2063SBarry Smith 
988cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
98917ab2063SBarry Smith {
990416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
991416022c9SBarry Smith   Scalar     *v = a->a;
99217ab2063SBarry Smith   double     sum = 0.0;
993416022c9SBarry Smith   int        i, j,shift = a->indexshift;
99417ab2063SBarry Smith 
99517ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
996416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
99717ab2063SBarry Smith #if defined(PETSC_COMPLEX)
99817ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
99917ab2063SBarry Smith #else
100017ab2063SBarry Smith       sum += (*v)*(*v); v++;
100117ab2063SBarry Smith #endif
100217ab2063SBarry Smith     }
100317ab2063SBarry Smith     *norm = sqrt(sum);
100417ab2063SBarry Smith   }
100517ab2063SBarry Smith   else if (type == NORM_1) {
100617ab2063SBarry Smith     double *tmp;
1007416022c9SBarry Smith     int    *jj = a->j;
10080452661fSBarry Smith     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
1009cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
101017ab2063SBarry Smith     *norm = 0.0;
1011416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1012a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
101317ab2063SBarry Smith     }
1014416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
101517ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
101617ab2063SBarry Smith     }
10170452661fSBarry Smith     PetscFree(tmp);
101817ab2063SBarry Smith   }
101917ab2063SBarry Smith   else if (type == NORM_INFINITY) {
102017ab2063SBarry Smith     *norm = 0.0;
1021416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1022416022c9SBarry Smith       v = a->a + a->i[j] + shift;
102317ab2063SBarry Smith       sum = 0.0;
1024416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1025cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
102617ab2063SBarry Smith       }
102717ab2063SBarry Smith       if (sum > *norm) *norm = sum;
102817ab2063SBarry Smith     }
102917ab2063SBarry Smith   }
103017ab2063SBarry Smith   else {
103148d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
103217ab2063SBarry Smith   }
103317ab2063SBarry Smith   return 0;
103417ab2063SBarry Smith }
103517ab2063SBarry Smith 
1036416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B)
103717ab2063SBarry Smith {
1038416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1039416022c9SBarry Smith   Mat        C;
1040416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1041416022c9SBarry Smith   int        shift = a->indexshift;
1042d5d45c9bSBarry Smith   Scalar     *array = a->a;
104317ab2063SBarry Smith 
10443638b69dSLois Curfman McInnes   if (B == PETSC_NULL && m != a->n)
1045dfa27b74SSatish Balay     SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place");
10460452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1047cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
104817ab2063SBarry Smith   if (shift) {
104917ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
105017ab2063SBarry Smith   }
105117ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1052416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
10530452661fSBarry Smith   PetscFree(col);
105417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
105517ab2063SBarry Smith     len = ai[i+1]-ai[i];
1056416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
105717ab2063SBarry Smith     array += len; aj += len;
105817ab2063SBarry Smith   }
105917ab2063SBarry Smith   if (shift) {
106017ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
106117ab2063SBarry Smith   }
106217ab2063SBarry Smith 
10636d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
10646d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
106517ab2063SBarry Smith 
10663638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1067416022c9SBarry Smith     *B = C;
106817ab2063SBarry Smith   } else {
1069416022c9SBarry Smith     /* This isn't really an in-place transpose */
10700452661fSBarry Smith     PetscFree(a->a);
10710452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
10720452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
10730452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
10740452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
10750452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
10761073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
10770452661fSBarry Smith     PetscFree(a);
1078416022c9SBarry Smith     PetscMemcpy(A,C,sizeof(struct _Mat));
10790452661fSBarry Smith     PetscHeaderDestroy(C);
108017ab2063SBarry Smith   }
108117ab2063SBarry Smith   return 0;
108217ab2063SBarry Smith }
108317ab2063SBarry Smith 
1084f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
108517ab2063SBarry Smith {
1086416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
108717ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1088d5d45c9bSBarry Smith   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
108917ab2063SBarry Smith 
109017ab2063SBarry Smith   if (ll) {
10913ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
10923ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
1093*9b1297e1SSatish Balay     VecGetLocalSize_Fast(ll,m);
1094f0b747eeSBarry Smith     if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length");
1095*9b1297e1SSatish Balay     VecGetArray_Fast(ll,l);
1096416022c9SBarry Smith     v = a->a;
109717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
109817ab2063SBarry Smith       x = l[i];
1099416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
110017ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
110117ab2063SBarry Smith     }
1102*9b1297e1SSatish Balay     VecRestoreArray_Fast(ll,l);
110344cd7ae7SLois Curfman McInnes     PLogFlops(nz);
110417ab2063SBarry Smith   }
110517ab2063SBarry Smith   if (rr) {
1106*9b1297e1SSatish Balay     VecGetLocalSize_Fast(rr,n);
1107f0b747eeSBarry Smith     if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length");
1108*9b1297e1SSatish Balay     VecGetArray_Fast(rr,r);
1109416022c9SBarry Smith     v = a->a; jj = a->j;
111017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
111117ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
111217ab2063SBarry Smith     }
1113*9b1297e1SSatish Balay     VecRestoreArray_Fast(rr,r);
111444cd7ae7SLois Curfman McInnes     PLogFlops(nz);
111517ab2063SBarry Smith   }
111617ab2063SBarry Smith   return 0;
111717ab2063SBarry Smith }
111817ab2063SBarry Smith 
1119cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
112017ab2063SBarry Smith {
1121db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
112202834360SBarry Smith   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
112399141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1124a2744918SBarry Smith   register int sum,lensi;
112599141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
112699141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
112799141d43SSatish Balay   Scalar       *a_new,*mat_a;
1128416022c9SBarry Smith   Mat          C;
112917ab2063SBarry Smith 
1130b48a1e75SSatish Balay   ierr = ISSorted(isrow,(PetscTruth*)&i);
1131b48a1e75SSatish Balay   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:ISrow is not sorted");
113299141d43SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
1133b48a1e75SSatish Balay   if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IScol is not sorted");
113499141d43SSatish Balay 
113517ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
113617ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
113717ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
113817ab2063SBarry Smith 
11397264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
114002834360SBarry Smith     /* special case of contiguous rows */
114157faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
114202834360SBarry Smith     starts = lens + ncols;
114302834360SBarry Smith     /* loop over new rows determining lens and starting points */
114402834360SBarry Smith     for (i=0; i<nrows; i++) {
1145a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1146a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
114702834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1148d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
114902834360SBarry Smith           starts[i] = k;
115002834360SBarry Smith           break;
115102834360SBarry Smith 	}
115202834360SBarry Smith       }
1153a2744918SBarry Smith       sum = 0;
115402834360SBarry Smith       while (k < kend) {
1155d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1156a2744918SBarry Smith         sum++;
115702834360SBarry Smith       }
1158a2744918SBarry Smith       lens[i] = sum;
115902834360SBarry Smith     }
116002834360SBarry Smith     /* create submatrix */
1161cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
116208480c60SBarry Smith       int n_cols,n_rows;
116308480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
116408480c60SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1165d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
116608480c60SBarry Smith       C = *B;
116708480c60SBarry Smith     }
116808480c60SBarry Smith     else {
116902834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
117008480c60SBarry Smith     }
1171db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1172db02288aSLois Curfman McInnes 
117302834360SBarry Smith     /* loop over rows inserting into submatrix */
1174db02288aSLois Curfman McInnes     a_new    = c->a;
1175db02288aSLois Curfman McInnes     j_new    = c->j;
1176db02288aSLois Curfman McInnes     i_new    = c->i;
1177db02288aSLois Curfman McInnes     i_new[0] = -shift;
117802834360SBarry Smith     for (i=0; i<nrows; i++) {
1179a2744918SBarry Smith       ii    = starts[i];
1180a2744918SBarry Smith       lensi = lens[i];
1181a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1182a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
118302834360SBarry Smith       }
1184a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1185a2744918SBarry Smith       a_new      += lensi;
1186a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1187a2744918SBarry Smith       c->ilen[i]  = lensi;
118802834360SBarry Smith     }
11890452661fSBarry Smith     PetscFree(lens);
119002834360SBarry Smith   }
119102834360SBarry Smith   else {
119202834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
11930452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
119402834360SBarry Smith     ssmap = smap + shift;
119599141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1196cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
119717ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
119802834360SBarry Smith     /* determine lens of each row */
119902834360SBarry Smith     for (i=0; i<nrows; i++) {
1200d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
120102834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
120202834360SBarry Smith       lens[i] = 0;
120302834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1204d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
120502834360SBarry Smith           lens[i]++;
120602834360SBarry Smith         }
120702834360SBarry Smith       }
120802834360SBarry Smith     }
120917ab2063SBarry Smith     /* Create and fill new matrix */
1210a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
121199141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
121299141d43SSatish Balay 
121399141d43SSatish Balay       if (c->m  != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
121499141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
121599141d43SSatish Balay         SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros");
121699141d43SSatish Balay       }
121799141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
121808480c60SBarry Smith       C = *B;
121908480c60SBarry Smith     }
122008480c60SBarry Smith     else {
122102834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
122208480c60SBarry Smith     }
122399141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
122417ab2063SBarry Smith     for (i=0; i<nrows; i++) {
122599141d43SSatish Balay       row    = irow[i];
122617ab2063SBarry Smith       nznew  = 0;
122799141d43SSatish Balay       kstart = ai[row]+shift;
122899141d43SSatish Balay       kend   = kstart + a->ilen[row];
122999141d43SSatish Balay       mat_i  = c->i[i]+shift;
123099141d43SSatish Balay       mat_j  = c->j + mat_i;
123199141d43SSatish Balay       mat_a  = c->a + mat_i;
123299141d43SSatish Balay       mat_ilen = c->ilen + i;
123317ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
123499141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
123599141d43SSatish Balay           *mat_j++ = tcol - (!shift);
123699141d43SSatish Balay           *mat_a++ = a->a[k];
123799141d43SSatish Balay           (*mat_ilen)++;
123899141d43SSatish Balay 
123917ab2063SBarry Smith         }
124017ab2063SBarry Smith       }
124117ab2063SBarry Smith     }
124202834360SBarry Smith     /* Free work space */
124302834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
124499141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
124502834360SBarry Smith   }
12466d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12476d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
124817ab2063SBarry Smith 
124917ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1250416022c9SBarry Smith   *B = C;
125117ab2063SBarry Smith   return 0;
125217ab2063SBarry Smith }
125317ab2063SBarry Smith 
1254a871dcd8SBarry Smith /*
125563b91edcSBarry Smith      note: This can only work for identity for row and col. It would
125663b91edcSBarry Smith    be good to check this and otherwise generate an error.
1257a871dcd8SBarry Smith */
125863b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1259a871dcd8SBarry Smith {
126063b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
126108480c60SBarry Smith   int        ierr;
126263b91edcSBarry Smith   Mat        outA;
126363b91edcSBarry Smith 
1264a871dcd8SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1265a871dcd8SBarry Smith 
126663b91edcSBarry Smith   outA          = inA;
126763b91edcSBarry Smith   inA->factor   = FACTOR_LU;
126863b91edcSBarry Smith   a->row        = row;
126963b91edcSBarry Smith   a->col        = col;
127063b91edcSBarry Smith 
12710452661fSBarry Smith   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
127263b91edcSBarry Smith 
127308480c60SBarry Smith   if (!a->diag) {
127408480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
127563b91edcSBarry Smith   }
127663b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1277a871dcd8SBarry Smith   return 0;
1278a871dcd8SBarry Smith }
1279a871dcd8SBarry Smith 
1280f0b747eeSBarry Smith #include "pinclude/plapack.h"
1281f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1282f0b747eeSBarry Smith {
1283f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1284f0b747eeSBarry Smith   int        one = 1;
1285f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1286f0b747eeSBarry Smith   PLogFlops(a->nz);
1287f0b747eeSBarry Smith   return 0;
1288f0b747eeSBarry Smith }
1289f0b747eeSBarry Smith 
1290cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1291cddf8d76SBarry Smith                                     Mat **B)
1292cddf8d76SBarry Smith {
1293cddf8d76SBarry Smith   int ierr,i;
1294cddf8d76SBarry Smith 
1295cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
12960452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1297cddf8d76SBarry Smith   }
1298cddf8d76SBarry Smith 
1299cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
1300905e6a2fSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
1301cddf8d76SBarry Smith   }
1302cddf8d76SBarry Smith   return 0;
1303cddf8d76SBarry Smith }
1304cddf8d76SBarry Smith 
13055a838052SSatish Balay static int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
13065a838052SSatish Balay {
13075a838052SSatish Balay   *bs = 1;
13085a838052SSatish Balay   return 0;
13095a838052SSatish Balay }
13105a838052SSatish Balay 
1311e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
13124dcbc457SBarry Smith {
1313e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
131406763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
13158a047759SSatish Balay   int        start, end, *ai, *aj;
131606763907SSatish Balay   char       *table;
13178a047759SSatish Balay   shift = a->indexshift;
1318e4d965acSSatish Balay   m     = a->m;
1319e4d965acSSatish Balay   ai    = a->i;
13208a047759SSatish Balay   aj    = a->j+shift;
13218a047759SSatish Balay 
1322e4d965acSSatish Balay   if (ov < 0)  SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used");
132306763907SSatish Balay 
132406763907SSatish Balay   table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table);
132506763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
132606763907SSatish Balay 
1327e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1328e4d965acSSatish Balay     /* Initialise the two local arrays */
1329e4d965acSSatish Balay     isz  = 0;
133006763907SSatish Balay     PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char));
1331e4d965acSSatish Balay 
1332e4d965acSSatish Balay                 /* Extract the indices, assume there can be duplicate entries */
13334dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
133477c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1335e4d965acSSatish Balay 
1336e4d965acSSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
1337e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
133806763907SSatish Balay       if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];}
13394dcbc457SBarry Smith     }
134006763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
134106763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1342e4d965acSSatish Balay 
134304a348a9SBarry Smith     k = 0;
134404a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap*/
134504a348a9SBarry Smith       n = isz;
134606763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1347e4d965acSSatish Balay         row   = nidx[k];
1348e4d965acSSatish Balay         start = ai[row];
1349e4d965acSSatish Balay         end   = ai[row+1];
135004a348a9SBarry Smith         for ( l = start; l<end ; l++){
13518a047759SSatish Balay           val = aj[l] + shift;
135206763907SSatish Balay           if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;}
1353e4d965acSSatish Balay         }
1354e4d965acSSatish Balay       }
1355e4d965acSSatish Balay     }
1356537820f0SBarry Smith     ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1357e4d965acSSatish Balay   }
135804a348a9SBarry Smith   PetscFree(table);
135906763907SSatish Balay   PetscFree(nidx);
1360e4d965acSSatish Balay   return 0;
13614dcbc457SBarry Smith }
136217ab2063SBarry Smith 
1363682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1364682d7d0cSBarry Smith {
1365682d7d0cSBarry Smith   static int called = 0;
1366682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1367682d7d0cSBarry Smith 
1368682d7d0cSBarry Smith   if (called) return 0; else called = 1;
136977c4ece6SBarry Smith   PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
137077c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_lu_pivotthreshold <threshold>\n");
137177c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n");
137277c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_no_inode  - Do not use inodes\n");
137377c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n");
1374682d7d0cSBarry Smith #if defined(HAVE_ESSL)
137577c4ece6SBarry Smith   PetscPrintf(comm,"  -mat_aij_essl  - Use IBM sparse LU factorization and solve.\n");
1376682d7d0cSBarry Smith #endif
1377682d7d0cSBarry Smith   return 0;
1378682d7d0cSBarry Smith }
1379205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1380a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1381a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1382a93ec695SBarry Smith 
1383682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
138417ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
138517ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1386416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1387416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
138817ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
138917ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
139017ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
139117ab2063SBarry Smith        MatRelax_SeqAIJ,
139217ab2063SBarry Smith        MatTranspose_SeqAIJ,
13937264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1394f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
139517ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
139617ab2063SBarry Smith        MatCompress_SeqAIJ,
139717ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
139817ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
139917ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
140017ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
140117ab2063SBarry Smith        0,0,MatConvert_SeqAIJ,
14023d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1403cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
14047eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1405682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1406f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
14075a838052SSatish Balay        MatScale_SeqAIJ,0,0,
14086945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
14096945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
14103b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
14113b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
14123b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1413a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1414a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
1415a93ec695SBarry Smith        MatColoringPatch_SeqAIJ};
141617ab2063SBarry Smith 
141717ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
141817ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
141917ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
142017ab2063SBarry Smith 
142117ab2063SBarry Smith /*@C
1422682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
14230d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
14246e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
14252bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
14262bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
142717ab2063SBarry Smith 
142817ab2063SBarry Smith    Input Parameters:
142917ab2063SBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
143017ab2063SBarry Smith .  m - number of rows
143117ab2063SBarry Smith .  n - number of columns
143217ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
14332bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
14342bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
143517ab2063SBarry Smith 
143617ab2063SBarry Smith    Output Parameter:
1437416022c9SBarry Smith .  A - the matrix
143817ab2063SBarry Smith 
143917ab2063SBarry Smith    Notes:
144017ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
144117ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
14420002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
144344cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
144417ab2063SBarry Smith 
144517ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1446a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
14473d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
14483d323bbdSBarry Smith    will get TERRIBLE performance, see the users' manual chapter on
14490d15e28bSLois Curfman McInnes    matrices and the file $(PETSC_DIR)/Performance.
145017ab2063SBarry Smith 
1451682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
1452682d7d0cSBarry Smith    improve numerical efficiency of Matrix vector products and solves. We
1453682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
14546c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
14556c7ebb05SLois Curfman McInnes 
14566c7ebb05SLois Curfman McInnes    Options Database Keys:
14576c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
14586e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
14596e62573dSLois Curfman McInnes $        (max limit=5)
14606e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
14616e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
14626e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
146317ab2063SBarry Smith 
146417ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
146517ab2063SBarry Smith @*/
1466416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
146717ab2063SBarry Smith {
1468416022c9SBarry Smith   Mat        B;
1469416022c9SBarry Smith   Mat_SeqAIJ *b;
14706945ee14SBarry Smith   int        i, len, ierr, flg,size;
14716945ee14SBarry Smith 
14726945ee14SBarry Smith   MPI_Comm_size(comm,&size);
14736945ee14SBarry Smith   if (size > 1) SETERRQ(1,"MatCreateSeqAIJ:Comm must be of size 1");
1474d5d45c9bSBarry Smith 
1475416022c9SBarry Smith   *A                  = 0;
14760452661fSBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1477416022c9SBarry Smith   PLogObjectCreate(B);
14780452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
147944cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1480416022c9SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1481416022c9SBarry Smith   B->destroy          = MatDestroy_SeqAIJ;
1482416022c9SBarry Smith   B->view             = MatView_SeqAIJ;
1483416022c9SBarry Smith   B->factor           = 0;
1484416022c9SBarry Smith   B->lupivotthreshold = 1.0;
14857a743949SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,
148669957df2SSatish Balay                           &flg); CHKERRQ(ierr);
14877a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
14887a743949SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",
14897a743949SBarry Smith                         (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr);
1490416022c9SBarry Smith   b->row              = 0;
1491416022c9SBarry Smith   b->col              = 0;
1492416022c9SBarry Smith   b->indexshift       = 0;
1493b810aeb4SBarry Smith   b->reallocs         = 0;
149469957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
149569957df2SSatish Balay   if (flg) b->indexshift = -1;
149617ab2063SBarry Smith 
149744cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
149844cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
14990452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1500b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
15017b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
15027b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1503416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
150417ab2063SBarry Smith     nz = nz*m;
150517ab2063SBarry Smith   }
150617ab2063SBarry Smith   else {
150717ab2063SBarry Smith     nz = 0;
1508416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
150917ab2063SBarry Smith   }
151017ab2063SBarry Smith 
151117ab2063SBarry Smith   /* allocate the matrix space */
1512416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
15130452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1514416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1515cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1516416022c9SBarry Smith   b->i  = b->j + nz;
1517416022c9SBarry Smith   b->singlemalloc = 1;
151817ab2063SBarry Smith 
1519416022c9SBarry Smith   b->i[0] = -b->indexshift;
152017ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1521416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
152217ab2063SBarry Smith   }
152317ab2063SBarry Smith 
1524416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
15250452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1526416022c9SBarry Smith   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1527416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
152817ab2063SBarry Smith 
1529416022c9SBarry Smith   b->nz               = 0;
1530416022c9SBarry Smith   b->maxnz            = nz;
1531416022c9SBarry Smith   b->sorted           = 0;
1532416022c9SBarry Smith   b->roworiented      = 1;
1533416022c9SBarry Smith   b->nonew            = 0;
1534416022c9SBarry Smith   b->diag             = 0;
1535416022c9SBarry Smith   b->solve_work       = 0;
153671bd300dSLois Curfman McInnes   b->spptr            = 0;
1537754ec7b1SSatish Balay   b->inode.node_count = 0;
1538754ec7b1SSatish Balay   b->inode.size       = 0;
15396c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
15406c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
15414e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
154217ab2063SBarry Smith 
1543416022c9SBarry Smith   *A = B;
15444e220ebcSLois Curfman McInnes 
15454b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
15464b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
154769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
154869957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
15494b14c69eSBarry Smith #endif
155069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
155169957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
155269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
155369957df2SSatish Balay   if (flg) {
1554416022c9SBarry Smith     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1555416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
155617ab2063SBarry Smith   }
155769957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
155869957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
155917ab2063SBarry Smith   return 0;
156017ab2063SBarry Smith }
156117ab2063SBarry Smith 
15623d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
156317ab2063SBarry Smith {
1564416022c9SBarry Smith   Mat        C;
1565416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
156608480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
156717ab2063SBarry Smith 
15684043dd9cSLois Curfman McInnes   *B = 0;
15690452661fSBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1570416022c9SBarry Smith   PLogObjectCreate(C);
15710452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
157241c01911SSatish Balay   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1573416022c9SBarry Smith   C->destroy    = MatDestroy_SeqAIJ;
1574416022c9SBarry Smith   C->view       = MatView_SeqAIJ;
1575416022c9SBarry Smith   C->factor     = A->factor;
1576416022c9SBarry Smith   c->row        = 0;
1577416022c9SBarry Smith   c->col        = 0;
1578416022c9SBarry Smith   c->indexshift = shift;
1579c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
158017ab2063SBarry Smith 
158144cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
158244cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
158344cd7ae7SLois Curfman McInnes   C->M          = a->m;
158444cd7ae7SLois Curfman McInnes   C->N          = a->n;
158517ab2063SBarry Smith 
15860452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
15870452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
158817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1589416022c9SBarry Smith     c->imax[i] = a->imax[i];
1590416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
159117ab2063SBarry Smith   }
159217ab2063SBarry Smith 
159317ab2063SBarry Smith   /* allocate the matrix space */
1594416022c9SBarry Smith   c->singlemalloc = 1;
1595416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
15960452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1597416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1598416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1599416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
160017ab2063SBarry Smith   if (m > 0) {
1601416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
160208480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1603416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
160417ab2063SBarry Smith     }
160508480c60SBarry Smith   }
160617ab2063SBarry Smith 
1607416022c9SBarry Smith   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1608416022c9SBarry Smith   c->sorted      = a->sorted;
1609416022c9SBarry Smith   c->roworiented = a->roworiented;
1610416022c9SBarry Smith   c->nonew       = a->nonew;
16117a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
161217ab2063SBarry Smith 
1613416022c9SBarry Smith   if (a->diag) {
16140452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1615416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
161617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1617416022c9SBarry Smith       c->diag[i] = a->diag[i];
161817ab2063SBarry Smith     }
161917ab2063SBarry Smith   }
1620416022c9SBarry Smith   else c->diag          = 0;
16216c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
16226c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1623754ec7b1SSatish Balay   if (a->inode.size){
1624daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1625754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1626daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1627754ec7b1SSatish Balay   } else {
1628754ec7b1SSatish Balay     c->inode.size       = 0;
1629754ec7b1SSatish Balay     c->inode.node_count = 0;
1630754ec7b1SSatish Balay   }
1631416022c9SBarry Smith   c->nz                 = a->nz;
1632416022c9SBarry Smith   c->maxnz              = a->maxnz;
1633416022c9SBarry Smith   c->solve_work         = 0;
163476dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1635754ec7b1SSatish Balay 
1636416022c9SBarry Smith   *B = C;
163717ab2063SBarry Smith   return 0;
163817ab2063SBarry Smith }
163917ab2063SBarry Smith 
164019bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
164117ab2063SBarry Smith {
1642416022c9SBarry Smith   Mat_SeqAIJ   *a;
1643416022c9SBarry Smith   Mat          B;
164417699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1645bcd2baecSBarry Smith   MPI_Comm     comm;
164617ab2063SBarry Smith 
164719bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
164817699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
164917699dbbSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
165090ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
165177c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
165248d91487SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
165317ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
165417ab2063SBarry Smith 
165517ab2063SBarry Smith   /* read in row lengths */
16560452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
165777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
165817ab2063SBarry Smith 
165917ab2063SBarry Smith   /* create our matrix */
1660416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1661416022c9SBarry Smith   B = *A;
1662416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1663416022c9SBarry Smith   shift = a->indexshift;
166417ab2063SBarry Smith 
166517ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
166677c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr);
166717ab2063SBarry Smith   if (shift) {
166817ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
1669416022c9SBarry Smith       a->j[i] += 1;
167017ab2063SBarry Smith     }
167117ab2063SBarry Smith   }
167217ab2063SBarry Smith 
167317ab2063SBarry Smith   /* read in nonzero values */
167477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr);
167517ab2063SBarry Smith 
167617ab2063SBarry Smith   /* set matrix "i" values */
1677416022c9SBarry Smith   a->i[0] = -shift;
167817ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
1679416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1680416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
168117ab2063SBarry Smith   }
16820452661fSBarry Smith   PetscFree(rowlengths);
168317ab2063SBarry Smith 
16846d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
16856d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
168617ab2063SBarry Smith   return 0;
168717ab2063SBarry Smith }
168817ab2063SBarry Smith 
1689686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
16907264ac53SSatish Balay {
16917264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
16927264ac53SSatish Balay 
1693bcd2baecSBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type");
16947264ac53SSatish Balay 
16957264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
16967264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
1697bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
169877c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1699bcd2baecSBarry Smith   }
17007264ac53SSatish Balay 
17017264ac53SSatish Balay   /* if the a->i are the same */
1702bcd2baecSBarry Smith   if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) {
170377c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
17047264ac53SSatish Balay   }
17057264ac53SSatish Balay 
17067264ac53SSatish Balay   /* if a->j are the same */
1707bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
170877c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
1709bcd2baecSBarry Smith   }
1710bcd2baecSBarry Smith 
1711bcd2baecSBarry Smith   /* if a->a are the same */
171219bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
171377c4ece6SBarry Smith     *flg = PETSC_FALSE; return 0;
17147264ac53SSatish Balay   }
171577c4ece6SBarry Smith   *flg = PETSC_TRUE;
17167264ac53SSatish Balay   return 0;
17177264ac53SSatish Balay 
17187264ac53SSatish Balay }
1719