xref: /petsc/src/mat/impls/aij/seq/aij.c (revision a5ae1ecdc7f32e56ab356dfc0ee0d859f7614a9d)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*a5ae1ecdSBarry Smith static char vcid[] = "$Id: aij.c,v 1.275 1998/07/14 02:35:46 bsmith Exp bsmith $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5d5d45c9bSBarry Smith /*
63369ce9aSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7d5d45c9bSBarry Smith   matrix storage format.
8d5d45c9bSBarry Smith */
93369ce9aSBarry Smith 
103369ce9aSBarry Smith #include "pinclude/pviewer.h"
113369ce9aSBarry Smith #include "sys.h"
1270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
13f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
14f5eb4b81SSatish Balay #include "src/inline/spops.h"
158d195f9aSBarry Smith #include "src/inline/dot.h"
16f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1717ab2063SBarry Smith 
18a2ce50c7SBarry Smith /*
19a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
20a2ce50c7SBarry Smith */
215615d1e5SSatish Balay #undef __FUNC__
225615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ"
23a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
24a2ce50c7SBarry Smith {
25a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
26a2ce50c7SBarry Smith   int        ierr = 1;
27a2ce50c7SBarry Smith 
283a40ed3dSBarry Smith   PetscFunctionBegin;
29e3372554SBarry Smith   SETERRQ(ierr,0,"Not implemented");
30d64ed03dSBarry Smith #if !defined(USE_PETSC_DEBUG)
31d64ed03dSBarry Smith   PetscFunctionReturn(0);
32d64ed03dSBarry Smith #endif
33a2ce50c7SBarry Smith }
34a2ce50c7SBarry Smith 
35bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3617ab2063SBarry Smith 
375615d1e5SSatish Balay #undef __FUNC__
385615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ"
398f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
406945ee14SBarry Smith                            PetscTruth *done)
4117ab2063SBarry Smith {
42416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
436945ee14SBarry Smith   int        ierr,i,ishift;
4417ab2063SBarry Smith 
453a40ed3dSBarry Smith   PetscFunctionBegin;
4631625ec5SSatish Balay   *m     = A->m;
473a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
486945ee14SBarry Smith   ishift = a->indexshift;
496945ee14SBarry Smith   if (symmetric) {
5031625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
516945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
5231625ec5SSatish Balay     int nz = a->i[a->m];
533b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
5431625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
553b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
563b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
5731625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
586945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
5931625ec5SSatish Balay     int nz = a->i[a->m] + 1;
603b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
6131625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
623b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
633b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
6431625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
656945ee14SBarry Smith   } else {
666945ee14SBarry Smith     *ia = a->i; *ja = a->j;
67a2ce50c7SBarry Smith   }
68a2ce50c7SBarry Smith 
693a40ed3dSBarry Smith   PetscFunctionReturn(0);
70a2744918SBarry Smith }
71a2744918SBarry Smith 
725615d1e5SSatish Balay #undef __FUNC__
735615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ"
748f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
756945ee14SBarry Smith                                PetscTruth *done)
766945ee14SBarry Smith {
776945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
783b2fbd54SBarry Smith   int        ishift = a->indexshift;
796945ee14SBarry Smith 
803a40ed3dSBarry Smith   PetscFunctionBegin;
813a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
823b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
836945ee14SBarry Smith     PetscFree(*ia);
846945ee14SBarry Smith     PetscFree(*ja);
85bcd2baecSBarry Smith   }
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
8717ab2063SBarry Smith }
8817ab2063SBarry Smith 
895615d1e5SSatish Balay #undef __FUNC__
905615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ"
9143a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
923b2fbd54SBarry Smith                            PetscTruth *done)
933b2fbd54SBarry Smith {
943b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
95a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
96a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
973b2fbd54SBarry Smith 
983a40ed3dSBarry Smith   PetscFunctionBegin;
993b2fbd54SBarry Smith   *nn     = A->n;
1003a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1013b2fbd54SBarry Smith   if (symmetric) {
102179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
1033b2fbd54SBarry Smith   } else {
10461d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
1053b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1063b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
107a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
1083b2fbd54SBarry Smith     jj = a->j;
1093b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
1103b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
1113b2fbd54SBarry Smith     }
1123b2fbd54SBarry Smith     cia[0] = oshift;
1133b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
1143b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
1153b2fbd54SBarry Smith     }
1163b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1173b2fbd54SBarry Smith     jj = a->j;
118a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
119a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
120a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1213b2fbd54SBarry Smith         col = *jj++ + ishift;
1223b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1233b2fbd54SBarry Smith       }
1243b2fbd54SBarry Smith     }
1253b2fbd54SBarry Smith     PetscFree(collengths);
1263b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1273b2fbd54SBarry Smith   }
1283b2fbd54SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1303b2fbd54SBarry Smith }
1313b2fbd54SBarry Smith 
1325615d1e5SSatish Balay #undef __FUNC__
1335615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ"
13443a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1353b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1363b2fbd54SBarry Smith {
1373a40ed3dSBarry Smith   PetscFunctionBegin;
1383a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1393b2fbd54SBarry Smith 
1403b2fbd54SBarry Smith   PetscFree(*ia);
1413b2fbd54SBarry Smith   PetscFree(*ja);
1423b2fbd54SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1443b2fbd54SBarry Smith }
1453b2fbd54SBarry Smith 
146227d817aSBarry Smith #define CHUNKSIZE   15
14717ab2063SBarry Smith 
1485615d1e5SSatish Balay #undef __FUNC__
1495615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ"
15061d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
15117ab2063SBarry Smith {
152416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
153416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1544b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
155d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
156416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
15717ab2063SBarry Smith 
1583a40ed3dSBarry Smith   PetscFunctionBegin;
15917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
160416022c9SBarry Smith     row  = im[k];
1613a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
162a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
163a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
1643b2fbd54SBarry Smith #endif
16517ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
16617ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
167416022c9SBarry Smith     low = 0;
16817ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
170a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
171a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
1723b2fbd54SBarry Smith #endif
1734b0e389bSBarry Smith       col = in[l] - shift;
1744b0e389bSBarry Smith       if (roworiented) {
1754b0e389bSBarry Smith         value = *v++;
176bef8e0ddSBarry Smith       } else {
1774b0e389bSBarry Smith         value = v[k + l*m];
1784b0e389bSBarry Smith       }
179416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
180416022c9SBarry Smith       while (high-low > 5) {
181416022c9SBarry Smith         t = (low+high)/2;
182416022c9SBarry Smith         if (rp[t] > col) high = t;
183416022c9SBarry Smith         else             low  = t;
18417ab2063SBarry Smith       }
185416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
18617ab2063SBarry Smith         if (rp[i] > col) break;
18717ab2063SBarry Smith         if (rp[i] == col) {
188416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
18917ab2063SBarry Smith           else                  ap[i] = value;
19017ab2063SBarry Smith           goto noinsert;
19117ab2063SBarry Smith         }
19217ab2063SBarry Smith       }
193c2653b3dSLois Curfman McInnes       if (nonew == 1) goto noinsert;
194a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
19517ab2063SBarry Smith       if (nrow >= rmax) {
19617ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
197416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
19817ab2063SBarry Smith         Scalar *new_a;
19917ab2063SBarry Smith 
200a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
20196854ed6SLois Curfman McInnes 
20217ab2063SBarry Smith         /* malloc new storage space */
203416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
2040452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
20517ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
20617ab2063SBarry Smith         new_i   = new_j + new_nz;
20717ab2063SBarry Smith 
20817ab2063SBarry Smith         /* copy over old data into new slots */
20917ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
210416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
211416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
212416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
213416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
21417ab2063SBarry Smith                                                            len*sizeof(int));
215416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
216416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
21717ab2063SBarry Smith                                                            len*sizeof(Scalar));
21817ab2063SBarry Smith         /* free up old matrix storage */
2190452661fSBarry Smith         PetscFree(a->a);
2200452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
221416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
222416022c9SBarry Smith         a->singlemalloc = 1;
22317ab2063SBarry Smith 
22417ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
225416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
226416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
227416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
228b810aeb4SBarry Smith         a->reallocs++;
22917ab2063SBarry Smith       }
230416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
231416022c9SBarry Smith       /* shift up all the later entries in this row */
232416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
23317ab2063SBarry Smith         rp[ii+1] = rp[ii];
23417ab2063SBarry Smith         ap[ii+1] = ap[ii];
23517ab2063SBarry Smith       }
23617ab2063SBarry Smith       rp[i] = col;
23717ab2063SBarry Smith       ap[i] = value;
23817ab2063SBarry Smith       noinsert:;
239416022c9SBarry Smith       low = i + 1;
24017ab2063SBarry Smith     }
24117ab2063SBarry Smith     ailen[row] = nrow;
24217ab2063SBarry Smith   }
2433a40ed3dSBarry Smith   PetscFunctionReturn(0);
24417ab2063SBarry Smith }
24517ab2063SBarry Smith 
2465615d1e5SSatish Balay #undef __FUNC__
2475615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ"
2488f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2497eb43aa7SLois Curfman McInnes {
2507eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
251b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2527eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2537eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2547eb43aa7SLois Curfman McInnes 
2553a40ed3dSBarry Smith   PetscFunctionBegin;
2567eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2577eb43aa7SLois Curfman McInnes     row  = im[k];
258a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
259a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
2607eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2617eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2627eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
263a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
264a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
2657eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2667eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2677eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2687eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2697eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2707eb43aa7SLois Curfman McInnes         else             low  = t;
2717eb43aa7SLois Curfman McInnes       }
2727eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2737eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2747eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
275b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2767eb43aa7SLois Curfman McInnes           goto finished;
2777eb43aa7SLois Curfman McInnes         }
2787eb43aa7SLois Curfman McInnes       }
279b49de8d1SLois Curfman McInnes       *v++ = zero;
2807eb43aa7SLois Curfman McInnes       finished:;
2817eb43aa7SLois Curfman McInnes     }
2827eb43aa7SLois Curfman McInnes   }
2833a40ed3dSBarry Smith   PetscFunctionReturn(0);
2847eb43aa7SLois Curfman McInnes }
2857eb43aa7SLois Curfman McInnes 
28617ab2063SBarry Smith 
2875615d1e5SSatish Balay #undef __FUNC__
2885615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary"
289480ef9eaSBarry Smith int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
29017ab2063SBarry Smith {
291416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
292416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
29317ab2063SBarry Smith 
2943a40ed3dSBarry Smith   PetscFunctionBegin;
29590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2960452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
297416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
298416022c9SBarry Smith   col_lens[1] = a->m;
299416022c9SBarry Smith   col_lens[2] = a->n;
300416022c9SBarry Smith   col_lens[3] = a->nz;
301416022c9SBarry Smith 
302416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
303416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
304416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
30517ab2063SBarry Smith   }
3060752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3070452661fSBarry Smith   PetscFree(col_lens);
308416022c9SBarry Smith 
309416022c9SBarry Smith   /* store column indices (zero start index) */
310416022c9SBarry Smith   if (a->indexshift) {
311416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
31217ab2063SBarry Smith   }
3130752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0); CHKERRQ(ierr);
314416022c9SBarry Smith   if (a->indexshift) {
315416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
31617ab2063SBarry Smith   }
317416022c9SBarry Smith 
318416022c9SBarry Smith   /* store nonzero values */
3190752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
32117ab2063SBarry Smith }
322416022c9SBarry Smith 
3235615d1e5SSatish Balay #undef __FUNC__
3245615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII"
325480ef9eaSBarry Smith int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
326416022c9SBarry Smith {
327416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
328496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
32917ab2063SBarry Smith   FILE        *fd;
33017ab2063SBarry Smith   char        *outputname;
33117ab2063SBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
33390ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
334416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
33590ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
336a93ec695SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
3373a40ed3dSBarry Smith     PetscFunctionReturn(0);
3383a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
339496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
340496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
341496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
34295e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
34395e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
3443a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
345d00d2cf4SBarry Smith     int nofinalvalue = 0;
346d00d2cf4SBarry Smith     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) {
347d00d2cf4SBarry Smith       nofinalvalue = 1;
348d00d2cf4SBarry Smith     }
349416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3504e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
351d00d2cf4SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue);
35217ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
35317ab2063SBarry Smith 
35417ab2063SBarry Smith     for (i=0; i<m; i++) {
355416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
3563a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
357e20fef11SSatish Balay         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
35817ab2063SBarry Smith #else
3597a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
36017ab2063SBarry Smith #endif
36117ab2063SBarry Smith       }
36217ab2063SBarry Smith     }
363d00d2cf4SBarry Smith     if (nofinalvalue) {
364d00d2cf4SBarry Smith       fprintf(fd,"%d %d  %18.16e\n", m, a->n, 0.0);
365d00d2cf4SBarry Smith     }
36617ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
3673a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
36844cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
36944cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
37044cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
3713a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
372e20fef11SSatish Balay         if (PetscImaginary(a->a[j]) > 0.0 && PetscReal(a->a[j]) != 0.0)
373e20fef11SSatish Balay           fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
374e20fef11SSatish Balay         else if (PetscImaginary(a->a[j]) < 0.0 && PetscReal(a->a[j]) != 0.0)
375e20fef11SSatish Balay           fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j]));
376e20fef11SSatish Balay         else if (PetscReal(a->a[j]) != 0.0)
377e20fef11SSatish Balay           fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j]));
37844cd7ae7SLois Curfman McInnes #else
37944cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
38044cd7ae7SLois Curfman McInnes #endif
38144cd7ae7SLois Curfman McInnes       }
38244cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
38344cd7ae7SLois Curfman McInnes     }
38444cd7ae7SLois Curfman McInnes   }
385496be53dSLois Curfman McInnes   else if (format == VIEWER_FORMAT_ASCII_SYMMODU) {
386496be53dSLois Curfman McInnes     int nzd=0, fshift=1, *sptr;
3872e44a96cSLois Curfman McInnes     sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr);
388496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
389496be53dSLois Curfman McInnes       sptr[i] = nzd+1;
390496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
391496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
3923a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
393e20fef11SSatish Balay           if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) nzd++;
394496be53dSLois Curfman McInnes #else
395496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) nzd++;
396496be53dSLois Curfman McInnes #endif
397496be53dSLois Curfman McInnes         }
398496be53dSLois Curfman McInnes       }
399496be53dSLois Curfman McInnes     }
4002e44a96cSLois Curfman McInnes     sptr[m] = nzd+1;
401496be53dSLois Curfman McInnes     fprintf(fd," %d %d\n\n",m,nzd);
4022e44a96cSLois Curfman McInnes     for ( i=0; i<m+1; i+=6 ) {
4032e44a96cSLois Curfman McInnes       if (i+4<m) fprintf(fd," %d %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);
4042e44a96cSLois Curfman McInnes       else if (i+3<m) fprintf(fd," %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);
4052e44a96cSLois Curfman McInnes       else if (i+2<m) fprintf(fd," %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);
4062e44a96cSLois Curfman McInnes       else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]);
4072e44a96cSLois Curfman McInnes       else if (i<m)   fprintf(fd," %d %d\n",sptr[i],sptr[i+1]);
4087272d637SLois Curfman McInnes       else            fprintf(fd," %d\n",sptr[i]);
409496be53dSLois Curfman McInnes     }
410496be53dSLois Curfman McInnes     fprintf(fd,"\n");
411496be53dSLois Curfman McInnes     PetscFree(sptr);
412496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
413496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
414496be53dSLois Curfman McInnes         if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift);
415496be53dSLois Curfman McInnes       }
416496be53dSLois Curfman McInnes       fprintf(fd,"\n");
417496be53dSLois Curfman McInnes     }
418496be53dSLois Curfman McInnes     fprintf(fd,"\n");
419496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
420496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
421496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
4223a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
423e20fef11SSatish Balay           if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0)
424e20fef11SSatish Balay             fprintf(fd," %18.16e %18.16e ",PetscReal(a->a[j]),PetscImaginary(a->a[j]));
425496be53dSLois Curfman McInnes #else
426496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]);
427496be53dSLois Curfman McInnes #endif
428496be53dSLois Curfman McInnes         }
429496be53dSLois Curfman McInnes       }
430496be53dSLois Curfman McInnes       fprintf(fd,"\n");
431496be53dSLois Curfman McInnes     }
4323a40ed3dSBarry Smith   } else {
43317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
43417ab2063SBarry Smith       fprintf(fd,"row %d:",i);
435416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
4363a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
437e20fef11SSatish Balay         if (PetscImaginary(a->a[j]) > 0.0) {
438e20fef11SSatish Balay           fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
439e20fef11SSatish Balay         } else if (PetscImaginary(a->a[j]) < 0.0) {
440e20fef11SSatish Balay           fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j]));
4413a40ed3dSBarry Smith         } else {
442e20fef11SSatish Balay           fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j]));
44317ab2063SBarry Smith         }
44417ab2063SBarry Smith #else
445416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
44617ab2063SBarry Smith #endif
44717ab2063SBarry Smith       }
44817ab2063SBarry Smith       fprintf(fd,"\n");
44917ab2063SBarry Smith     }
45017ab2063SBarry Smith   }
45117ab2063SBarry Smith   fflush(fd);
4523a40ed3dSBarry Smith   PetscFunctionReturn(0);
453416022c9SBarry Smith }
454416022c9SBarry Smith 
4555615d1e5SSatish Balay #undef __FUNC__
456480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw_Zoom"
457480ef9eaSBarry Smith int MatView_SeqAIJ_Draw_Zoom(Draw draw,void *Aa)
458416022c9SBarry Smith {
459480ef9eaSBarry Smith   Mat         A = (Mat) Aa;
460416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
461480ef9eaSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,color;
4620513a670SBarry Smith   int         format;
463480ef9eaSBarry Smith   double      xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
464480ef9eaSBarry Smith   Viewer      viewer;
465cddf8d76SBarry Smith 
4663a40ed3dSBarry Smith   PetscFunctionBegin;
467480ef9eaSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr);
4680513a670SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
46919bcc07fSBarry Smith 
470480ef9eaSBarry Smith   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr);
471416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
4720513a670SBarry Smith 
4730513a670SBarry Smith   if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4740513a670SBarry Smith     /* Blue for negative, Cyan for zero and  Red for positive */
475cddf8d76SBarry Smith     color = DRAW_BLUE;
476416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
477cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
478416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
479cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
4803a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
481e20fef11SSatish Balay         if (PetscReal(a->a[j]) >=  0.) continue;
482cddf8d76SBarry Smith #else
483cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
484cddf8d76SBarry Smith #endif
485480ef9eaSBarry Smith         ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
486cddf8d76SBarry Smith       }
487cddf8d76SBarry Smith     }
488cddf8d76SBarry Smith     color = DRAW_CYAN;
489cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
490cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
491cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
492cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
493cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
494480ef9eaSBarry Smith         ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
495cddf8d76SBarry Smith       }
496cddf8d76SBarry Smith     }
497cddf8d76SBarry Smith     color = DRAW_RED;
498cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
499cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
500cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
501cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5023a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
503e20fef11SSatish Balay         if (PetscReal(a->a[j]) <=  0.) continue;
504cddf8d76SBarry Smith #else
505cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
506cddf8d76SBarry Smith #endif
507480ef9eaSBarry Smith         ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
508416022c9SBarry Smith       }
509416022c9SBarry Smith     }
5100513a670SBarry Smith   } else {
5110513a670SBarry Smith     /* use contour shading to indicate magnitude of values */
5120513a670SBarry Smith     /* first determine max of all nonzero values */
5130513a670SBarry Smith     int    nz = a->nz,count;
5140513a670SBarry Smith     Draw   popup;
515480ef9eaSBarry Smith     double scale;
5160513a670SBarry Smith 
5170513a670SBarry Smith     for ( i=0; i<nz; i++ ) {
5180513a670SBarry Smith       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
5190513a670SBarry Smith     }
520480ef9eaSBarry Smith     scale = (245.0 - DRAW_BASIC_COLORS)/maxv;
521522c5e43SBarry Smith     ierr  = DrawGetPopup(draw,&popup); CHKERRQ(ierr);
5220513a670SBarry Smith     ierr  = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr);
5230513a670SBarry Smith     count = 0;
5240513a670SBarry Smith     for ( i=0; i<m; i++ ) {
5250513a670SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
5260513a670SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5270513a670SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5286d6b6c30SSatish Balay         color = DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(a->a[count]));
529480ef9eaSBarry Smith         ierr  = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
5300513a670SBarry Smith         count++;
5310513a670SBarry Smith       }
5320513a670SBarry Smith     }
5330513a670SBarry Smith   }
534480ef9eaSBarry Smith   PetscFunctionReturn(0);
535480ef9eaSBarry Smith }
536cddf8d76SBarry Smith 
537480ef9eaSBarry Smith #undef __FUNC__
538480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw"
539480ef9eaSBarry Smith int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
540480ef9eaSBarry Smith {
541480ef9eaSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data;
542480ef9eaSBarry Smith   int        ierr;
543480ef9eaSBarry Smith   Draw       draw;
544480ef9eaSBarry Smith   double     xr,yr,xl,yl,h,w;
545480ef9eaSBarry Smith 
546480ef9eaSBarry Smith   PetscTruth isnull;
547480ef9eaSBarry Smith 
548480ef9eaSBarry Smith   PetscFunctionBegin;
549480ef9eaSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
550480ef9eaSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr);
551480ef9eaSBarry Smith   if (isnull) PetscFunctionReturn(0);
552480ef9eaSBarry Smith 
553480ef9eaSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
554480ef9eaSBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
555480ef9eaSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
556cddf8d76SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
557480ef9eaSBarry Smith   ierr = DrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A); CHKERRQ(ierr);
558480ef9eaSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5593a40ed3dSBarry Smith   PetscFunctionReturn(0);
560416022c9SBarry Smith }
561416022c9SBarry Smith 
5625615d1e5SSatish Balay #undef __FUNC__
563d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ"
564e1311b90SBarry Smith int MatView_SeqAIJ(Mat A,Viewer viewer)
565416022c9SBarry Smith {
566416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
567bcd2baecSBarry Smith   ViewerType  vtype;
568bcd2baecSBarry Smith   int         ierr;
569416022c9SBarry Smith 
5703a40ed3dSBarry Smith   PetscFunctionBegin;
571bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
572bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
5733a40ed3dSBarry Smith     ierr = ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);  CHKERRQ(ierr);
5745cd90555SBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
5753a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_ASCII(A,viewer); CHKERRQ(ierr);
5765cd90555SBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
5773a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_Binary(A,viewer); CHKERRQ(ierr);
5785cd90555SBarry Smith   } else if (vtype == DRAW_VIEWER) {
5793a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_Draw(A,viewer); CHKERRQ(ierr);
5805cd90555SBarry Smith   } else {
5815cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
58217ab2063SBarry Smith   }
5833a40ed3dSBarry Smith   PetscFunctionReturn(0);
58417ab2063SBarry Smith }
58519bcc07fSBarry Smith 
586c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
5875615d1e5SSatish Balay #undef __FUNC__
5885615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ"
5898f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
59017ab2063SBarry Smith {
591416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
59241c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
59343ee02c3SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0;
594416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
59517ab2063SBarry Smith 
5963a40ed3dSBarry Smith   PetscFunctionBegin;
5973a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
59817ab2063SBarry Smith 
59943ee02c3SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
60017ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
601416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
60217ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
60394a9d846SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
60417ab2063SBarry Smith     if (fshift) {
605416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
60617ab2063SBarry Smith       N = ailen[i];
60717ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
60817ab2063SBarry Smith         ip[j-fshift] = ip[j];
60917ab2063SBarry Smith         ap[j-fshift] = ap[j];
61017ab2063SBarry Smith       }
61117ab2063SBarry Smith     }
61217ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
61317ab2063SBarry Smith   }
61417ab2063SBarry Smith   if (m) {
61517ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
61617ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
61717ab2063SBarry Smith   }
61817ab2063SBarry Smith   /* reset ilen and imax for each row */
61917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
62017ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
62117ab2063SBarry Smith   }
622416022c9SBarry Smith   a->nz = ai[m] + shift;
62317ab2063SBarry Smith 
62417ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
625416022c9SBarry Smith   if (fshift && a->diag) {
6260452661fSBarry Smith     PetscFree(a->diag);
627416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
628416022c9SBarry Smith     a->diag = 0;
62917ab2063SBarry Smith   }
6304e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
6314e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
6324e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
633b810aeb4SBarry Smith            a->reallocs);
63494a9d846SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax);
635dd5f02e7SSatish Balay   a->reallocs          = 0;
6364e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
6374e220ebcSLois Curfman McInnes 
63876dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
63941c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
6403a40ed3dSBarry Smith   PetscFunctionReturn(0);
64117ab2063SBarry Smith }
64217ab2063SBarry Smith 
6435615d1e5SSatish Balay #undef __FUNC__
6445615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ"
6458f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A)
64617ab2063SBarry Smith {
647416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
6483a40ed3dSBarry Smith 
6493a40ed3dSBarry Smith   PetscFunctionBegin;
650cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
6513a40ed3dSBarry Smith   PetscFunctionReturn(0);
65217ab2063SBarry Smith }
653416022c9SBarry Smith 
6545615d1e5SSatish Balay #undef __FUNC__
6555615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ"
656e1311b90SBarry Smith int MatDestroy_SeqAIJ(Mat A)
65717ab2063SBarry Smith {
658416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
65982bf6240SBarry Smith   int        ierr;
660d5d45c9bSBarry Smith 
6613a40ed3dSBarry Smith   PetscFunctionBegin;
66294d884c6SBarry Smith   if (--A->refct > 0) PetscFunctionReturn(0);
66394d884c6SBarry Smith 
66494d884c6SBarry Smith   if (A->mapping) {
66594d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr);
66694d884c6SBarry Smith   }
66794d884c6SBarry Smith   if (A->bmapping) {
66894d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr);
66994d884c6SBarry Smith   }
67061b13de0SBarry Smith   if (A->rmap) {
67161b13de0SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
67261b13de0SBarry Smith   }
67361b13de0SBarry Smith   if (A->cmap) {
67461b13de0SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
67561b13de0SBarry Smith   }
6763a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
677e1311b90SBarry Smith   PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
67817ab2063SBarry Smith #endif
6790452661fSBarry Smith   PetscFree(a->a);
6800452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
6810452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
6820452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
6830452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
6840452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
68576dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
68682bf6240SBarry Smith   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
6870452661fSBarry Smith   PetscFree(a);
688eed86810SBarry Smith 
689f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
690f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
6913a40ed3dSBarry Smith   PetscFunctionReturn(0);
69217ab2063SBarry Smith }
69317ab2063SBarry Smith 
6945615d1e5SSatish Balay #undef __FUNC__
6955615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
6968f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A)
69717ab2063SBarry Smith {
6983a40ed3dSBarry Smith   PetscFunctionBegin;
6993a40ed3dSBarry Smith   PetscFunctionReturn(0);
70017ab2063SBarry Smith }
70117ab2063SBarry Smith 
7025615d1e5SSatish Balay #undef __FUNC__
7035615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
7048f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op)
70517ab2063SBarry Smith {
706416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7073a40ed3dSBarry Smith 
7083a40ed3dSBarry Smith   PetscFunctionBegin;
7096d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
7106d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
7116d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
712219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
7136d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
714c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
71596854ed6SLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
7166d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
7176d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
718219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7196d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7206d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
72190f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
722b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES||
723b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
724981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n");
7253a40ed3dSBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS) {
7263a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7273a40ed3dSBarry Smith   } else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
7286d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
7296d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
7306d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
7316d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
7323a40ed3dSBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7333a40ed3dSBarry Smith   PetscFunctionReturn(0);
73417ab2063SBarry Smith }
73517ab2063SBarry Smith 
7365615d1e5SSatish Balay #undef __FUNC__
7375615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
7388f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
73917ab2063SBarry Smith {
740416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7413a40ed3dSBarry Smith   int        i,j, n,shift = a->indexshift,ierr;
74217ab2063SBarry Smith   Scalar     *x, zero = 0.0;
74317ab2063SBarry Smith 
7443a40ed3dSBarry Smith   PetscFunctionBegin;
7453a40ed3dSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
746e1311b90SBarry Smith   ierr = VecGetArray(v,&x);;CHKERRQ(ierr);
747e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr);
748a8c6a408SBarry Smith   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
749416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
750416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
751416022c9SBarry Smith       if (a->j[j]+shift == i) {
752416022c9SBarry Smith         x[i] = a->a[j];
75317ab2063SBarry Smith         break;
75417ab2063SBarry Smith       }
75517ab2063SBarry Smith     }
75617ab2063SBarry Smith   }
757e1311b90SBarry Smith   ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr);
7583a40ed3dSBarry Smith   PetscFunctionReturn(0);
75917ab2063SBarry Smith }
76017ab2063SBarry Smith 
76117ab2063SBarry Smith /* -------------------------------------------------------*/
76217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
76317ab2063SBarry Smith /* -------------------------------------------------------*/
7645615d1e5SSatish Balay #undef __FUNC__
7655615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
76644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
76717ab2063SBarry Smith {
768416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
76917ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
770e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx, shift = a->indexshift;
77117ab2063SBarry Smith 
7723a40ed3dSBarry Smith   PetscFunctionBegin;
773e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
774e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
775cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
77617ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
77717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
778416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
779416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
780416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
78117ab2063SBarry Smith     alpha = x[i];
78217ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
78317ab2063SBarry Smith   }
784416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
785e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
786e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7873a40ed3dSBarry Smith   PetscFunctionReturn(0);
78817ab2063SBarry Smith }
789d5d45c9bSBarry Smith 
7905615d1e5SSatish Balay #undef __FUNC__
7915615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
79244cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
79317ab2063SBarry Smith {
794416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
79517ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
796e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx,shift = a->indexshift;
79717ab2063SBarry Smith 
7983a40ed3dSBarry Smith   PetscFunctionBegin;
799e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
800e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
80117ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
80217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
80317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
804416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
805416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
806416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
80717ab2063SBarry Smith     alpha = x[i];
80817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
80917ab2063SBarry Smith   }
81090f02eecSBarry Smith   PLogFlops(2*a->nz);
811e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
812e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8133a40ed3dSBarry Smith   PetscFunctionReturn(0);
81417ab2063SBarry Smith }
81517ab2063SBarry Smith 
8165615d1e5SSatish Balay #undef __FUNC__
8175615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
81844cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
81917ab2063SBarry Smith {
820416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
82117ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
822e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
8232e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTAIJ)
824e36a17ebSSatish Balay   int        n, i, jrow,j;
825e36a17ebSSatish Balay #endif
82617ab2063SBarry Smith 
827fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT)
828fee21e36SBarry Smith #pragma disjoint(*x,*y,*v)
829fee21e36SBarry Smith #endif
830fee21e36SBarry Smith 
8313a40ed3dSBarry Smith   PetscFunctionBegin;
832e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
833e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
83417ab2063SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
835416022c9SBarry Smith   idx  = a->j;
836d64ed03dSBarry Smith   v    = a->a;
837416022c9SBarry Smith   ii   = a->i;
838acc4d009SSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTAIJ)
83929b5ca8aSSatish Balay   fortranmultaij_(&m,x,ii,idx+shift,v+shift,y);
8408d195f9aSBarry Smith #else
841d64ed03dSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
842d64ed03dSBarry Smith   idx  += shift;
84317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8449ea0dfa2SSatish Balay     jrow = ii[i];
8459ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
84617ab2063SBarry Smith     sum  = 0.0;
8479ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8489ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8499ea0dfa2SSatish Balay      }
85017ab2063SBarry Smith     y[i] = sum;
85117ab2063SBarry Smith   }
8528d195f9aSBarry Smith #endif
853416022c9SBarry Smith   PLogFlops(2*a->nz - m);
854e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
855e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
8563a40ed3dSBarry Smith   PetscFunctionReturn(0);
85717ab2063SBarry Smith }
85817ab2063SBarry Smith 
8595615d1e5SSatish Balay #undef __FUNC__
8605615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
86144cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
86217ab2063SBarry Smith {
863416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
86417ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
865e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
8662e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTADDAIJ)
867e36a17ebSSatish Balay   int        n,i,jrow,j;
868e36a17ebSSatish Balay #endif
8699ea0dfa2SSatish Balay 
8703a40ed3dSBarry Smith   PetscFunctionBegin;
871e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
872e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
873e1311b90SBarry Smith   ierr = VecGetArray(zz,&z); CHKERRQ(ierr);
87417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
875cddf8d76SBarry Smith   idx  = a->j;
876d64ed03dSBarry Smith   v    = a->a;
877cddf8d76SBarry Smith   ii   = a->i;
8782e40c62eSSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTADDAIJ)
87929b5ca8aSSatish Balay   fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
88002ab625aSSatish Balay #else
881d64ed03dSBarry Smith   v   += shift; /* shift for Fortran start by 1 indexing */
882d64ed03dSBarry Smith   idx += shift;
88317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8849ea0dfa2SSatish Balay     jrow = ii[i];
8859ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
88617ab2063SBarry Smith     sum  = y[i];
8879ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8889ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8899ea0dfa2SSatish Balay      }
89017ab2063SBarry Smith     z[i] = sum;
89117ab2063SBarry Smith   }
89202ab625aSSatish Balay #endif
893416022c9SBarry Smith   PLogFlops(2*a->nz);
894e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
895e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
896e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr);
8973a40ed3dSBarry Smith   PetscFunctionReturn(0);
89817ab2063SBarry Smith }
89917ab2063SBarry Smith 
90017ab2063SBarry Smith /*
90117ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
90217ab2063SBarry Smith */
90317ab2063SBarry Smith 
9045615d1e5SSatish Balay #undef __FUNC__
9055615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
906416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
90717ab2063SBarry Smith {
908416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
909416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
91017ab2063SBarry Smith 
9113a40ed3dSBarry Smith   PetscFunctionBegin;
9120452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
913416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
914416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
915416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
916416022c9SBarry Smith       if (a->j[j]+shift == i) {
91717ab2063SBarry Smith         diag[i] = j - shift;
91817ab2063SBarry Smith         break;
91917ab2063SBarry Smith       }
92017ab2063SBarry Smith     }
92117ab2063SBarry Smith   }
922416022c9SBarry Smith   a->diag = diag;
9233a40ed3dSBarry Smith   PetscFunctionReturn(0);
92417ab2063SBarry Smith }
92517ab2063SBarry Smith 
9265615d1e5SSatish Balay #undef __FUNC__
9275615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
92844cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
92917ab2063SBarry Smith                            double fshift,int its,Vec xx)
93017ab2063SBarry Smith {
931416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
932416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
933d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
93417ab2063SBarry Smith 
9353a40ed3dSBarry Smith   PetscFunctionBegin;
936e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
937e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
938d00d2cf4SBarry Smith   if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);}
939416022c9SBarry Smith   diag = a->diag;
940416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
94117ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
94217ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
94317ab2063SBarry Smith     bs = b + shift;
94417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
945416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
946416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
947416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
948416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
94917ab2063SBarry Smith         sum  = b[i]*d/omega;
95017ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
95117ab2063SBarry Smith         x[i] = sum;
95217ab2063SBarry Smith     }
9533a40ed3dSBarry Smith     PetscFunctionReturn(0);
95417ab2063SBarry Smith   }
95517ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
956a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done");
9573a40ed3dSBarry Smith   } else if (flag & SOR_EISENSTAT) {
95817ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
95917ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
96017ab2063SBarry Smith 
96117ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
96217ab2063SBarry Smith 
96317ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
96417ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
96517ab2063SBarry Smith     is the relaxation factor.
96617ab2063SBarry Smith     */
9670452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
96817ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
96917ab2063SBarry Smith 
97017ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
97117ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
972416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
973416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
974416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
975416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
97617ab2063SBarry Smith       sum  = b[i];
97717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
97817ab2063SBarry Smith       x[i] = omega*(sum/d);
97917ab2063SBarry Smith     }
98017ab2063SBarry Smith 
98117ab2063SBarry Smith     /*  t = b - (2*E - D)x */
982416022c9SBarry Smith     v = a->a;
98317ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
98417ab2063SBarry Smith 
98517ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
986416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
987416022c9SBarry Smith     diag = a->diag;
98817ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
989416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
990416022c9SBarry Smith       n    = diag[i] - a->i[i];
991416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
992416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
99317ab2063SBarry Smith       sum  = t[i];
99417ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
99517ab2063SBarry Smith       t[i] = omega*(sum/d);
99617ab2063SBarry Smith     }
99717ab2063SBarry Smith 
99817ab2063SBarry Smith     /*  x = x + t */
99917ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
10000452661fSBarry Smith     PetscFree(t);
10013a40ed3dSBarry Smith     PetscFunctionReturn(0);
100217ab2063SBarry Smith   }
100317ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
100417ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
100517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1006416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1007416022c9SBarry Smith         n    = diag[i] - a->i[i];
1008416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1009416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
101017ab2063SBarry Smith         sum  = b[i];
101117ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
101217ab2063SBarry Smith         x[i] = omega*(sum/d);
101317ab2063SBarry Smith       }
101417ab2063SBarry Smith       xb = x;
10153a40ed3dSBarry Smith     } else xb = b;
101617ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
101717ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
101817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1019416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
102017ab2063SBarry Smith       }
102117ab2063SBarry Smith     }
102217ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
102317ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1024416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1025416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
1026416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1027416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
102817ab2063SBarry Smith         sum  = xb[i];
102917ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
103017ab2063SBarry Smith         x[i] = omega*(sum/d);
103117ab2063SBarry Smith       }
103217ab2063SBarry Smith     }
103317ab2063SBarry Smith     its--;
103417ab2063SBarry Smith   }
103517ab2063SBarry Smith   while (its--) {
103617ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
103717ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1038416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1039416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1040416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1041416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
104217ab2063SBarry Smith         sum  = b[i];
104317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10447e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
104517ab2063SBarry Smith       }
104617ab2063SBarry Smith     }
104717ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
104817ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1049416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1050416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1051416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1052416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
105317ab2063SBarry Smith         sum  = b[i];
105417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10557e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
105617ab2063SBarry Smith       }
105717ab2063SBarry Smith     }
105817ab2063SBarry Smith   }
10593a40ed3dSBarry Smith   PetscFunctionReturn(0);
106017ab2063SBarry Smith }
106117ab2063SBarry Smith 
10625615d1e5SSatish Balay #undef __FUNC__
10635615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10648f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
106517ab2063SBarry Smith {
1066416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
10674e220ebcSLois Curfman McInnes 
10683a40ed3dSBarry Smith   PetscFunctionBegin;
10694e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
10704e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
10714e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
10724e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
10734e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10744e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
10754e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
10764e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
10774e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
10784e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
10794e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
10804e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
10814e220ebcSLois Curfman McInnes   info->memory         = A->mem;
10824e220ebcSLois Curfman McInnes   if (A->factor) {
10834e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10844e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
10854e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
10864e220ebcSLois Curfman McInnes   } else {
10874e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
10884e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
10894e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
10904e220ebcSLois Curfman McInnes   }
10913a40ed3dSBarry Smith   PetscFunctionReturn(0);
109217ab2063SBarry Smith }
109317ab2063SBarry Smith 
109417ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
109517ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
109617ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
109717ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
109817ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
109917ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
110017ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
110117ab2063SBarry Smith 
11025615d1e5SSatish Balay #undef __FUNC__
11035615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
11048f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
110517ab2063SBarry Smith {
1106416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1107416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
110817ab2063SBarry Smith 
11093a40ed3dSBarry Smith   PetscFunctionBegin;
111077c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
111117ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
111217ab2063SBarry Smith   if (diag) {
111317ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1114a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1115416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1116416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1117416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1118416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
11193a40ed3dSBarry Smith       } else {
1120d64ed03dSBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr);
112117ab2063SBarry Smith       }
112217ab2063SBarry Smith     }
11233a40ed3dSBarry Smith   } else {
112417ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1125a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1126416022c9SBarry Smith       a->ilen[rows[i]] = 0;
112717ab2063SBarry Smith     }
112817ab2063SBarry Smith   }
112917ab2063SBarry Smith   ISRestoreIndices(is,&rows);
113043a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11313a40ed3dSBarry Smith   PetscFunctionReturn(0);
113217ab2063SBarry Smith }
113317ab2063SBarry Smith 
11345615d1e5SSatish Balay #undef __FUNC__
11355615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11368f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
113717ab2063SBarry Smith {
1138416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11393a40ed3dSBarry Smith 
11403a40ed3dSBarry Smith   PetscFunctionBegin;
11410752156aSBarry Smith   if (m) *m = a->m;
11420752156aSBarry Smith   if (n) *n = a->n;
11433a40ed3dSBarry Smith   PetscFunctionReturn(0);
114417ab2063SBarry Smith }
114517ab2063SBarry Smith 
11465615d1e5SSatish Balay #undef __FUNC__
11475615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11488f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
114917ab2063SBarry Smith {
1150416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11513a40ed3dSBarry Smith 
11523a40ed3dSBarry Smith   PetscFunctionBegin;
1153416022c9SBarry Smith   *m = 0; *n = a->m;
11543a40ed3dSBarry Smith   PetscFunctionReturn(0);
115517ab2063SBarry Smith }
1156026e39d0SSatish Balay 
11575615d1e5SSatish Balay #undef __FUNC__
11585615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11594e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
116017ab2063SBarry Smith {
1161416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1162c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
116317ab2063SBarry Smith 
11643a40ed3dSBarry Smith   PetscFunctionBegin;
1165a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
116617ab2063SBarry Smith 
1167416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1168416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
116917ab2063SBarry Smith   if (idx) {
1170416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
11714e093b46SBarry Smith     if (*nz && shift) {
11720452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
117317ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
11744e093b46SBarry Smith     } else if (*nz) {
11754e093b46SBarry Smith       *idx = itmp;
117617ab2063SBarry Smith     }
117717ab2063SBarry Smith     else *idx = 0;
117817ab2063SBarry Smith   }
11793a40ed3dSBarry Smith   PetscFunctionReturn(0);
118017ab2063SBarry Smith }
118117ab2063SBarry Smith 
11825615d1e5SSatish Balay #undef __FUNC__
11835615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
11844e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
118517ab2063SBarry Smith {
11864e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11873a40ed3dSBarry Smith 
11883a40ed3dSBarry Smith   PetscFunctionBegin;
11894e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
11903a40ed3dSBarry Smith   PetscFunctionReturn(0);
119117ab2063SBarry Smith }
119217ab2063SBarry Smith 
11935615d1e5SSatish Balay #undef __FUNC__
11945615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
11958f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
119617ab2063SBarry Smith {
1197416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1198416022c9SBarry Smith   Scalar     *v = a->a;
119917ab2063SBarry Smith   double     sum = 0.0;
1200416022c9SBarry Smith   int        i, j,shift = a->indexshift;
120117ab2063SBarry Smith 
12023a40ed3dSBarry Smith   PetscFunctionBegin;
120317ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1204416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
12053a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
1206e20fef11SSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
120717ab2063SBarry Smith #else
120817ab2063SBarry Smith       sum += (*v)*(*v); v++;
120917ab2063SBarry Smith #endif
121017ab2063SBarry Smith     }
121117ab2063SBarry Smith     *norm = sqrt(sum);
12123a40ed3dSBarry Smith   } else if (type == NORM_1) {
121317ab2063SBarry Smith     double *tmp;
1214416022c9SBarry Smith     int    *jj = a->j;
121566963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1216cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
121717ab2063SBarry Smith     *norm = 0.0;
1218416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1219a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
122017ab2063SBarry Smith     }
1221416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
122217ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
122317ab2063SBarry Smith     }
12240452661fSBarry Smith     PetscFree(tmp);
12253a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
122617ab2063SBarry Smith     *norm = 0.0;
1227416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1228416022c9SBarry Smith       v = a->a + a->i[j] + shift;
122917ab2063SBarry Smith       sum = 0.0;
1230416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1231cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
123217ab2063SBarry Smith       }
123317ab2063SBarry Smith       if (sum > *norm) *norm = sum;
123417ab2063SBarry Smith     }
12353a40ed3dSBarry Smith   } else {
1236a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
123717ab2063SBarry Smith   }
12383a40ed3dSBarry Smith   PetscFunctionReturn(0);
123917ab2063SBarry Smith }
124017ab2063SBarry Smith 
12415615d1e5SSatish Balay #undef __FUNC__
12425615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12438f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
124417ab2063SBarry Smith {
1245416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1246416022c9SBarry Smith   Mat        C;
1247416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1248416022c9SBarry Smith   int        shift = a->indexshift;
1249d5d45c9bSBarry Smith   Scalar     *array = a->a;
125017ab2063SBarry Smith 
12513a40ed3dSBarry Smith   PetscFunctionBegin;
1252a8c6a408SBarry Smith   if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place");
12530452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1254cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
125517ab2063SBarry Smith   if (shift) {
125617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
125717ab2063SBarry Smith   }
125817ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1259416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12600452661fSBarry Smith   PetscFree(col);
126117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
126217ab2063SBarry Smith     len = ai[i+1]-ai[i];
1263416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
126417ab2063SBarry Smith     array += len; aj += len;
126517ab2063SBarry Smith   }
126617ab2063SBarry Smith   if (shift) {
126717ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
126817ab2063SBarry Smith   }
126917ab2063SBarry Smith 
12706d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12716d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
127217ab2063SBarry Smith 
12733638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1274416022c9SBarry Smith     *B = C;
127517ab2063SBarry Smith   } else {
1276f830108cSBarry Smith     PetscOps *Abops;
12770a6ffc59SBarry Smith     MatOps   Aops;
1278f830108cSBarry Smith 
1279416022c9SBarry Smith     /* This isn't really an in-place transpose */
12800452661fSBarry Smith     PetscFree(a->a);
12810452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
12820452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
12830452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
12840452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
12850452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
12861073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
12870452661fSBarry Smith     PetscFree(a);
1288f830108cSBarry Smith 
1289f830108cSBarry Smith     /*
1290f830108cSBarry Smith       This is horrible, horrible code. We need to keep the
12918f0f457cSSatish Balay       the bops and ops Structures, copy everything from C
12928f0f457cSSatish Balay       including the function pointers..
1293f830108cSBarry Smith     */
12948f0f457cSSatish Balay     PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps));
12958f0f457cSSatish Balay     PetscMemcpy(A->bops,C->bops,sizeof(PetscOps));
1296f830108cSBarry Smith     Abops    = A->bops;
1297f830108cSBarry Smith     Aops     = A->ops;
1298f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
1299f830108cSBarry Smith     A->bops  = Abops;
1300f830108cSBarry Smith     A->ops   = Aops;
130127a8da17SBarry Smith     A->qlist = 0;
1302f830108cSBarry Smith 
13030452661fSBarry Smith     PetscHeaderDestroy(C);
130417ab2063SBarry Smith   }
13053a40ed3dSBarry Smith   PetscFunctionReturn(0);
130617ab2063SBarry Smith }
130717ab2063SBarry Smith 
13085615d1e5SSatish Balay #undef __FUNC__
13095615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
13108f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
131117ab2063SBarry Smith {
1312416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
131317ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1314e1311b90SBarry Smith   int        ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
131517ab2063SBarry Smith 
13163a40ed3dSBarry Smith   PetscFunctionBegin;
131717ab2063SBarry Smith   if (ll) {
13183ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
13193ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
1320e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr);
1321a8c6a408SBarry Smith     if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
1322e1311b90SBarry Smith     ierr = VecGetArray(ll,&l); CHKERRQ(ierr);
1323416022c9SBarry Smith     v = a->a;
132417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
132517ab2063SBarry Smith       x = l[i];
1326416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
132717ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
132817ab2063SBarry Smith     }
1329e1311b90SBarry Smith     ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr);
133044cd7ae7SLois Curfman McInnes     PLogFlops(nz);
133117ab2063SBarry Smith   }
133217ab2063SBarry Smith   if (rr) {
1333e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr);
1334a8c6a408SBarry Smith     if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
1335e1311b90SBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1336416022c9SBarry Smith     v = a->a; jj = a->j;
133717ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
133817ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
133917ab2063SBarry Smith     }
1340e1311b90SBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
134144cd7ae7SLois Curfman McInnes     PLogFlops(nz);
134217ab2063SBarry Smith   }
13433a40ed3dSBarry Smith   PetscFunctionReturn(0);
134417ab2063SBarry Smith }
134517ab2063SBarry Smith 
13465615d1e5SSatish Balay #undef __FUNC__
13475615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
13486a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B)
134917ab2063SBarry Smith {
1350db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1351d5db1b72SBarry Smith   int          *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
135299141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1353a2744918SBarry Smith   register int sum,lensi;
135499141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
135599141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
135699141d43SSatish Balay   Scalar       *a_new,*mat_a;
1357416022c9SBarry Smith   Mat          C;
1358fee21e36SBarry Smith   PetscTruth   stride;
135917ab2063SBarry Smith 
13603a40ed3dSBarry Smith   PetscFunctionBegin;
1361d64ed03dSBarry Smith   ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr);
1362a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
1363d64ed03dSBarry Smith   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
1364a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
136599141d43SSatish Balay 
136617ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
136717ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
136817ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
136917ab2063SBarry Smith 
1370fee21e36SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr);
1371fee21e36SBarry Smith   ierr = ISStride(iscol,&stride); CHKERRQ(ierr);
1372fee21e36SBarry Smith   if (stride && step == 1) {
137302834360SBarry Smith     /* special case of contiguous rows */
137457faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
137502834360SBarry Smith     starts = lens + ncols;
137602834360SBarry Smith     /* loop over new rows determining lens and starting points */
137702834360SBarry Smith     for (i=0; i<nrows; i++) {
1378a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1379a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
138002834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1381d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
138202834360SBarry Smith           starts[i] = k;
138302834360SBarry Smith           break;
138402834360SBarry Smith 	}
138502834360SBarry Smith       }
1386a2744918SBarry Smith       sum = 0;
138702834360SBarry Smith       while (k < kend) {
1388d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1389a2744918SBarry Smith         sum++;
139002834360SBarry Smith       }
1391a2744918SBarry Smith       lens[i] = sum;
139202834360SBarry Smith     }
139302834360SBarry Smith     /* create submatrix */
1394cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
139508480c60SBarry Smith       int n_cols,n_rows;
139608480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1397a8c6a408SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1398d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
139908480c60SBarry Smith       C = *B;
14003a40ed3dSBarry Smith     } else {
140102834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
140208480c60SBarry Smith     }
1403db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1404db02288aSLois Curfman McInnes 
140502834360SBarry Smith     /* loop over rows inserting into submatrix */
1406db02288aSLois Curfman McInnes     a_new    = c->a;
1407db02288aSLois Curfman McInnes     j_new    = c->j;
1408db02288aSLois Curfman McInnes     i_new    = c->i;
1409db02288aSLois Curfman McInnes     i_new[0] = -shift;
141002834360SBarry Smith     for (i=0; i<nrows; i++) {
1411a2744918SBarry Smith       ii    = starts[i];
1412a2744918SBarry Smith       lensi = lens[i];
1413a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1414a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
141502834360SBarry Smith       }
1416a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1417a2744918SBarry Smith       a_new      += lensi;
1418a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1419a2744918SBarry Smith       c->ilen[i]  = lensi;
142002834360SBarry Smith     }
14210452661fSBarry Smith     PetscFree(lens);
14223a40ed3dSBarry Smith   } else {
142302834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
14240452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
142502834360SBarry Smith     ssmap = smap + shift;
142699141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1427cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
142817ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
142902834360SBarry Smith     /* determine lens of each row */
143002834360SBarry Smith     for (i=0; i<nrows; i++) {
1431d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
143202834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
143302834360SBarry Smith       lens[i] = 0;
143402834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1435d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
143602834360SBarry Smith           lens[i]++;
143702834360SBarry Smith         }
143802834360SBarry Smith       }
143902834360SBarry Smith     }
144017ab2063SBarry Smith     /* Create and fill new matrix */
1441a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
144299141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
1443a8c6a408SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
144499141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1445a8c6a408SBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
144699141d43SSatish Balay       }
144799141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
144808480c60SBarry Smith       C = *B;
14493a40ed3dSBarry Smith     } else {
145002834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
145108480c60SBarry Smith     }
145299141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
145317ab2063SBarry Smith     for (i=0; i<nrows; i++) {
145499141d43SSatish Balay       row    = irow[i];
145599141d43SSatish Balay       kstart = ai[row]+shift;
145699141d43SSatish Balay       kend   = kstart + a->ilen[row];
145799141d43SSatish Balay       mat_i  = c->i[i]+shift;
145899141d43SSatish Balay       mat_j  = c->j + mat_i;
145999141d43SSatish Balay       mat_a  = c->a + mat_i;
146099141d43SSatish Balay       mat_ilen = c->ilen + i;
146117ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
146299141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
146399141d43SSatish Balay           *mat_j++ = tcol - (!shift);
146499141d43SSatish Balay           *mat_a++ = a->a[k];
146599141d43SSatish Balay           (*mat_ilen)++;
146699141d43SSatish Balay 
146717ab2063SBarry Smith         }
146817ab2063SBarry Smith       }
146917ab2063SBarry Smith     }
147002834360SBarry Smith     /* Free work space */
147102834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
147299141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
147302834360SBarry Smith   }
14746d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
14756d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
147617ab2063SBarry Smith 
147717ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1478416022c9SBarry Smith   *B = C;
14793a40ed3dSBarry Smith   PetscFunctionReturn(0);
148017ab2063SBarry Smith }
148117ab2063SBarry Smith 
1482a871dcd8SBarry Smith /*
148363b91edcSBarry Smith      note: This can only work for identity for row and col. It would
148463b91edcSBarry Smith    be good to check this and otherwise generate an error.
1485a871dcd8SBarry Smith */
14865615d1e5SSatish Balay #undef __FUNC__
14875615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
14888f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1489a871dcd8SBarry Smith {
149063b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
149108480c60SBarry Smith   int        ierr;
149263b91edcSBarry Smith   Mat        outA;
149363b91edcSBarry Smith 
14943a40ed3dSBarry Smith   PetscFunctionBegin;
1495a8c6a408SBarry Smith   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
1496a871dcd8SBarry Smith 
149763b91edcSBarry Smith   outA          = inA;
149863b91edcSBarry Smith   inA->factor   = FACTOR_LU;
149963b91edcSBarry Smith   a->row        = row;
150063b91edcSBarry Smith   a->col        = col;
150163b91edcSBarry Smith 
1502f0ec6fceSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1503f0ec6fceSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
15041e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1505f0ec6fceSSatish Balay 
150694a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
15070452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work);
150894a9d846SBarry Smith   }
150963b91edcSBarry Smith 
151008480c60SBarry Smith   if (!a->diag) {
151108480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
151263b91edcSBarry Smith   }
151363b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
15143a40ed3dSBarry Smith   PetscFunctionReturn(0);
1515a871dcd8SBarry Smith }
1516a871dcd8SBarry Smith 
151774cdf0dfSBarry Smith #include "pinclude/blaslapack.h"
15185615d1e5SSatish Balay #undef __FUNC__
15195615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
15208f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1521f0b747eeSBarry Smith {
1522f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1523f0b747eeSBarry Smith   int        one = 1;
15243a40ed3dSBarry Smith 
15253a40ed3dSBarry Smith   PetscFunctionBegin;
1526f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1527f0b747eeSBarry Smith   PLogFlops(a->nz);
15283a40ed3dSBarry Smith   PetscFunctionReturn(0);
1529f0b747eeSBarry Smith }
1530f0b747eeSBarry Smith 
15315615d1e5SSatish Balay #undef __FUNC__
15325615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
15336a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1534cddf8d76SBarry Smith {
1535cddf8d76SBarry Smith   int ierr,i;
1536cddf8d76SBarry Smith 
15373a40ed3dSBarry Smith   PetscFunctionBegin;
1538cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
15390452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1540cddf8d76SBarry Smith   }
1541cddf8d76SBarry Smith 
1542cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
15436a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1544cddf8d76SBarry Smith   }
15453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1546cddf8d76SBarry Smith }
1547cddf8d76SBarry Smith 
15485615d1e5SSatish Balay #undef __FUNC__
15495615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
15508f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
15515a838052SSatish Balay {
1552f830108cSBarry Smith   PetscFunctionBegin;
15535a838052SSatish Balay   *bs = 1;
15543a40ed3dSBarry Smith   PetscFunctionReturn(0);
15555a838052SSatish Balay }
15565a838052SSatish Balay 
15575615d1e5SSatish Balay #undef __FUNC__
15585615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15598f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15604dcbc457SBarry Smith {
1561e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
156206763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15638a047759SSatish Balay   int        start, end, *ai, *aj;
1564bbd702dbSSatish Balay   BT         table;
1565bbd702dbSSatish Balay 
15663a40ed3dSBarry Smith   PetscFunctionBegin;
15678a047759SSatish Balay   shift = a->indexshift;
1568e4d965acSSatish Balay   m     = a->m;
1569e4d965acSSatish Balay   ai    = a->i;
15708a047759SSatish Balay   aj    = a->j+shift;
15718a047759SSatish Balay 
1572a8c6a408SBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used");
157306763907SSatish Balay 
157406763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1575bbd702dbSSatish Balay   ierr  = BTCreate(m,table); CHKERRQ(ierr);
157606763907SSatish Balay 
1577e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1578b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1579e4d965acSSatish Balay     isz  = 0;
1580bbd702dbSSatish Balay     BTMemzero(m,table);
1581e4d965acSSatish Balay 
1582e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
15834dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
158477c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1585e4d965acSSatish Balay 
1586dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1587e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
1588bbd702dbSSatish Balay       if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];}
15894dcbc457SBarry Smith     }
159006763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
159106763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1592e4d965acSSatish Balay 
159304a348a9SBarry Smith     k = 0;
159404a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
159504a348a9SBarry Smith       n = isz;
159606763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1597e4d965acSSatish Balay         row   = nidx[k];
1598e4d965acSSatish Balay         start = ai[row];
1599e4d965acSSatish Balay         end   = ai[row+1];
160004a348a9SBarry Smith         for ( l = start; l<end ; l++){
16018a047759SSatish Balay           val = aj[l] + shift;
16022a8f2162SSatish Balay           if (!BTLookupSet(table,val)) {nidx[isz++] = val;}
1603e4d965acSSatish Balay         }
1604e4d965acSSatish Balay       }
1605e4d965acSSatish Balay     }
1606029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1607e4d965acSSatish Balay   }
1608bbd702dbSSatish Balay   BTDestroy(table);
160906763907SSatish Balay   PetscFree(nidx);
16103a40ed3dSBarry Smith   PetscFunctionReturn(0);
16114dcbc457SBarry Smith }
161217ab2063SBarry Smith 
16130513a670SBarry Smith /* -------------------------------------------------------------- */
16140513a670SBarry Smith #undef __FUNC__
16150513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
16160513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
16170513a670SBarry Smith {
16180513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
16190513a670SBarry Smith   Scalar     *vwork;
16200513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
16210513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
162256cd22aeSBarry Smith   IS         icolp,irowp;
16230513a670SBarry Smith 
16243a40ed3dSBarry Smith   PetscFunctionBegin;
162556cd22aeSBarry Smith   ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr);
162656cd22aeSBarry Smith   ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr);
162756cd22aeSBarry Smith   ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr);
162856cd22aeSBarry Smith   ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr);
16290513a670SBarry Smith 
16300513a670SBarry Smith   /* determine lengths of permuted rows */
16310513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
16320513a670SBarry Smith   for (i=0; i<m; i++ ) {
16330513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
16340513a670SBarry Smith   }
16350513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
16360513a670SBarry Smith   PetscFree(lens);
16370513a670SBarry Smith 
16380513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
16390513a670SBarry Smith   for (i=0; i<m; i++) {
16400513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16410513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
16420513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
16430513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16440513a670SBarry Smith   }
16450513a670SBarry Smith   PetscFree(cnew);
16460513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
16470513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
164856cd22aeSBarry Smith   ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr);
164956cd22aeSBarry Smith   ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr);
165056cd22aeSBarry Smith   ierr = ISDestroy(irowp); CHKERRQ(ierr);
165156cd22aeSBarry Smith   ierr = ISDestroy(icolp); CHKERRQ(ierr);
16523a40ed3dSBarry Smith   PetscFunctionReturn(0);
16530513a670SBarry Smith }
16540513a670SBarry Smith 
16555615d1e5SSatish Balay #undef __FUNC__
16565615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1657682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1658682d7d0cSBarry Smith {
1659682d7d0cSBarry Smith   static int called = 0;
1660682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1661682d7d0cSBarry Smith 
16623a40ed3dSBarry Smith   PetscFunctionBegin;
16633a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
166476be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
166576be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
166676be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
166776be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodes\n");
166876be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1669682d7d0cSBarry Smith #if defined(HAVE_ESSL)
167076be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1671682d7d0cSBarry Smith #endif
16723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1673682d7d0cSBarry Smith }
16748f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1675a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1676a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1677a93ec695SBarry Smith 
1678682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
16790a6ffc59SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
168017ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1681416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1682416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
168317ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
168417ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
168517ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
168617ab2063SBarry Smith        MatRelax_SeqAIJ,
168717ab2063SBarry Smith        MatTranspose_SeqAIJ,
16887264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1689f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
169017ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
169117ab2063SBarry Smith        MatCompress_SeqAIJ,
169217ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
169317ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
169417ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
169517ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
169694a9d846SBarry Smith        0,0,
16973d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1698cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
16997eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1700682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1701f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
17025a838052SSatish Balay        MatScale_SeqAIJ,0,0,
17036945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
17046945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
17053b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
17063b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
17073b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1708a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1709a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
17100513a670SBarry Smith        MatColoringPatch_SeqAIJ,
17110513a670SBarry Smith        0,
1712cda55fadSBarry Smith        MatPermute_SeqAIJ,
1713cda55fadSBarry Smith        0,
1714cda55fadSBarry Smith        0,
1715cda55fadSBarry Smith        0,
1716cda55fadSBarry Smith        0,
1717cda55fadSBarry Smith        MatGetMaps_Petsc};
171817ab2063SBarry Smith 
171917ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
172017ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
172117ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
172217ab2063SBarry Smith 
17235615d1e5SSatish Balay #undef __FUNC__
1724bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ"
1725bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
1726bef8e0ddSBarry Smith {
1727bef8e0ddSBarry Smith   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
1728bef8e0ddSBarry Smith   int        i,nz,n;
1729bef8e0ddSBarry Smith 
1730bef8e0ddSBarry Smith   PetscFunctionBegin;
1731bef8e0ddSBarry Smith   if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing");
1732bef8e0ddSBarry Smith 
1733bef8e0ddSBarry Smith   nz = aij->maxnz;
1734bef8e0ddSBarry Smith   n  = aij->n;
1735bef8e0ddSBarry Smith   for (i=0; i<nz; i++) {
1736bef8e0ddSBarry Smith     aij->j[i] = indices[i];
1737bef8e0ddSBarry Smith   }
1738bef8e0ddSBarry Smith   aij->nz = nz;
1739bef8e0ddSBarry Smith   for ( i=0; i<n; i++ ) {
1740bef8e0ddSBarry Smith     aij->ilen[i] = aij->imax[i];
1741bef8e0ddSBarry Smith   }
1742bef8e0ddSBarry Smith 
1743bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1744bef8e0ddSBarry Smith }
1745bef8e0ddSBarry Smith 
1746bef8e0ddSBarry Smith #undef __FUNC__
1747bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices"
1748bef8e0ddSBarry Smith /*@
1749bef8e0ddSBarry Smith     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
1750bef8e0ddSBarry Smith        in the matrix.
1751bef8e0ddSBarry Smith 
1752bef8e0ddSBarry Smith   Input Parameters:
1753bef8e0ddSBarry Smith +  mat - the SeqAIJ matrix
1754bef8e0ddSBarry Smith -  indices - the column indices
1755bef8e0ddSBarry Smith 
1756bef8e0ddSBarry Smith   Notes:
1757bef8e0ddSBarry Smith     This can be called if you have precomputed the nonzero structure of the
1758bef8e0ddSBarry Smith   matrix and want to provide it to the matrix object to improve the performance
1759bef8e0ddSBarry Smith   of the MatSetValues() operation.
1760bef8e0ddSBarry Smith 
1761bef8e0ddSBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
1762bef8e0ddSBarry Smith   MatCreateSeqAIJ().
1763bef8e0ddSBarry Smith 
1764bef8e0ddSBarry Smith     MUST be called before any calls to MatSetValues();
1765bef8e0ddSBarry Smith 
1766bef8e0ddSBarry Smith @*/
1767bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
1768bef8e0ddSBarry Smith {
1769bef8e0ddSBarry Smith   int ierr,(*f)(Mat,int *);
1770bef8e0ddSBarry Smith 
1771bef8e0ddSBarry Smith   PetscFunctionBegin;
1772bef8e0ddSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1773bef8e0ddSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f);
1774bef8e0ddSBarry Smith          CHKERRQ(ierr);
1775bef8e0ddSBarry Smith   if (f) {
1776bef8e0ddSBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
1777bef8e0ddSBarry Smith   } else {
1778bef8e0ddSBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
1779bef8e0ddSBarry Smith   }
1780bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1781bef8e0ddSBarry Smith }
1782bef8e0ddSBarry Smith 
1783bef8e0ddSBarry Smith #undef __FUNC__
17845615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
178517ab2063SBarry Smith /*@C
1786682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
17870d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
17886e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
17892bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
17902bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
179117ab2063SBarry Smith 
1792db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1793db81eaa0SLois Curfman McInnes 
179417ab2063SBarry Smith    Input Parameters:
1795db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
179617ab2063SBarry Smith .  m - number of rows
179717ab2063SBarry Smith .  n - number of columns
179817ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
1799db81eaa0SLois Curfman McInnes -  nzz - array containing the number of nonzeros in the various rows
18002bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
180117ab2063SBarry Smith 
180217ab2063SBarry Smith    Output Parameter:
1803416022c9SBarry Smith .  A - the matrix
180417ab2063SBarry Smith 
1805b259b22eSLois Curfman McInnes    Notes:
180617ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
180717ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
18080002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
180944cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
181017ab2063SBarry Smith 
181117ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1812a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
18133d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
18146da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
181517ab2063SBarry Smith 
1816682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
18174fca80b9SLois Curfman McInnes    improve numerical efficiency of matrix-vector products and solves. We
1818682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
18196c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
18206c7ebb05SLois Curfman McInnes 
18216c7ebb05SLois Curfman McInnes    Options Database Keys:
1822db81eaa0SLois Curfman McInnes +  -mat_aij_no_inode  - Do not use inodes
1823db81eaa0SLois Curfman McInnes .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
1824db81eaa0SLois Curfman McInnes -  -mat_aij_oneindex - Internally use indexing starting at 1
1825db81eaa0SLois Curfman McInnes         rather than 0.  Note that when calling MatSetValues(),
1826db81eaa0SLois Curfman McInnes         the user still MUST index entries starting at 0!
182717ab2063SBarry Smith 
1828bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices()
182917ab2063SBarry Smith @*/
1830416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
183117ab2063SBarry Smith {
1832416022c9SBarry Smith   Mat        B;
1833416022c9SBarry Smith   Mat_SeqAIJ *b;
18346945ee14SBarry Smith   int        i, len, ierr, flg,size;
18356945ee14SBarry Smith 
18363a40ed3dSBarry Smith   PetscFunctionBegin;
18376945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1838a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1");
1839d5d45c9bSBarry Smith 
1840416022c9SBarry Smith   *A                  = 0;
1841f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView);
1842416022c9SBarry Smith   PLogObjectCreate(B);
18430452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
184444cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
18450a6ffc59SBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1846e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_SeqAIJ;
1847e1311b90SBarry Smith   B->ops->view             = MatView_SeqAIJ;
1848416022c9SBarry Smith   B->factor           = 0;
1849416022c9SBarry Smith   B->lupivotthreshold = 1.0;
185090f02eecSBarry Smith   B->mapping          = 0;
1851e1311b90SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg); CHKERRQ(ierr);
18527a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
1853e1311b90SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*) &b->ilu_preserve_row_sums);CHKERRQ(ierr);
1854416022c9SBarry Smith   b->row              = 0;
1855416022c9SBarry Smith   b->col              = 0;
185682bf6240SBarry Smith   b->icol             = 0;
1857416022c9SBarry Smith   b->indexshift       = 0;
1858b810aeb4SBarry Smith   b->reallocs         = 0;
185969957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
186069957df2SSatish Balay   if (flg) b->indexshift = -1;
186117ab2063SBarry Smith 
186244cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
186344cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
1864*a5ae1ecdSBarry Smith 
1865*a5ae1ecdSBarry Smith   ierr = MapCreate(comm,m,m,B->rmap);CHKERRQ(ierr);
1866*a5ae1ecdSBarry Smith   ierr = MapCreate(comm,n,n,B->cmap);CHKERRQ(ierr);
1867*a5ae1ecdSBarry Smith 
18680452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1869b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
18707b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
18717b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1872416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
187317ab2063SBarry Smith     nz = nz*m;
18743a40ed3dSBarry Smith   } else {
187517ab2063SBarry Smith     nz = 0;
1876416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
187717ab2063SBarry Smith   }
187817ab2063SBarry Smith 
187917ab2063SBarry Smith   /* allocate the matrix space */
1880416022c9SBarry Smith   len             = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
18810452661fSBarry Smith   b->a            = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1882416022c9SBarry Smith   b->j            = (int *) (b->a + nz);
1883cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1884416022c9SBarry Smith   b->i            = b->j + nz;
1885416022c9SBarry Smith   b->singlemalloc = 1;
188617ab2063SBarry Smith 
1887416022c9SBarry Smith   b->i[0] = -b->indexshift;
188817ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1889416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
189017ab2063SBarry Smith   }
189117ab2063SBarry Smith 
1892416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
18930452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1894f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1895416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
189617ab2063SBarry Smith 
1897416022c9SBarry Smith   b->nz               = 0;
1898416022c9SBarry Smith   b->maxnz            = nz;
1899416022c9SBarry Smith   b->sorted           = 0;
1900416022c9SBarry Smith   b->roworiented      = 1;
1901416022c9SBarry Smith   b->nonew            = 0;
1902416022c9SBarry Smith   b->diag             = 0;
1903416022c9SBarry Smith   b->solve_work       = 0;
190471bd300dSLois Curfman McInnes   b->spptr            = 0;
1905754ec7b1SSatish Balay   b->inode.node_count = 0;
1906754ec7b1SSatish Balay   b->inode.size       = 0;
19076c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
19086c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
19094e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
191017ab2063SBarry Smith 
1911416022c9SBarry Smith   *A = B;
19124e220ebcSLois Curfman McInnes 
19134b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
19144b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
191569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
191669957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
19174b14c69eSBarry Smith #endif
191869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
191969957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
192069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
192169957df2SSatish Balay   if (flg) {
1922a8c6a408SBarry Smith     if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1923416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
192417ab2063SBarry Smith   }
192569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
192669957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1927bef8e0ddSBarry Smith 
1928bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
1929bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
1930bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
19313a40ed3dSBarry Smith   PetscFunctionReturn(0);
193217ab2063SBarry Smith }
193317ab2063SBarry Smith 
19345615d1e5SSatish Balay #undef __FUNC__
19355615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
19363d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
193717ab2063SBarry Smith {
1938416022c9SBarry Smith   Mat        C;
1939416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1940bef8e0ddSBarry Smith   int        i,len, m = a->m,shift = a->indexshift,ierr;
194117ab2063SBarry Smith 
19423a40ed3dSBarry Smith   PetscFunctionBegin;
19434043dd9cSLois Curfman McInnes   *B = 0;
1944f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView);
1945416022c9SBarry Smith   PLogObjectCreate(C);
19460452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1947f830108cSBarry Smith   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1948e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqAIJ;
1949e1311b90SBarry Smith   C->ops->view       = MatView_SeqAIJ;
1950416022c9SBarry Smith   C->factor     = A->factor;
1951416022c9SBarry Smith   c->row        = 0;
1952416022c9SBarry Smith   c->col        = 0;
195382bf6240SBarry Smith   c->icol       = 0;
1954416022c9SBarry Smith   c->indexshift = shift;
1955c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
195617ab2063SBarry Smith 
195744cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
195844cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
195944cd7ae7SLois Curfman McInnes   C->M          = a->m;
196044cd7ae7SLois Curfman McInnes   C->N          = a->n;
196117ab2063SBarry Smith 
19620452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
19630452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
196417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1965416022c9SBarry Smith     c->imax[i] = a->imax[i];
1966416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
196717ab2063SBarry Smith   }
196817ab2063SBarry Smith 
196917ab2063SBarry Smith   /* allocate the matrix space */
1970416022c9SBarry Smith   c->singlemalloc = 1;
1971416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
19720452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1973416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1974416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1975416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
197617ab2063SBarry Smith   if (m > 0) {
1977416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
197808480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1979416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
198017ab2063SBarry Smith     }
198108480c60SBarry Smith   }
198217ab2063SBarry Smith 
1983f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1984416022c9SBarry Smith   c->sorted      = a->sorted;
1985416022c9SBarry Smith   c->roworiented = a->roworiented;
1986416022c9SBarry Smith   c->nonew       = a->nonew;
19877a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
198817ab2063SBarry Smith 
1989416022c9SBarry Smith   if (a->diag) {
19900452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1991416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
199217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1993416022c9SBarry Smith       c->diag[i] = a->diag[i];
199417ab2063SBarry Smith     }
19953a40ed3dSBarry Smith   } else c->diag          = 0;
19966c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
19976c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1998754ec7b1SSatish Balay   if (a->inode.size){
1999daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
2000754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
2001daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
2002754ec7b1SSatish Balay   } else {
2003754ec7b1SSatish Balay     c->inode.size       = 0;
2004754ec7b1SSatish Balay     c->inode.node_count = 0;
2005754ec7b1SSatish Balay   }
2006416022c9SBarry Smith   c->nz                 = a->nz;
2007416022c9SBarry Smith   c->maxnz              = a->maxnz;
2008416022c9SBarry Smith   c->solve_work         = 0;
200976dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
2010754ec7b1SSatish Balay 
2011416022c9SBarry Smith   *B = C;
2012bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqAIJSetColumnIndices_C",
2013bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2014bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
20153a40ed3dSBarry Smith   PetscFunctionReturn(0);
201617ab2063SBarry Smith }
201717ab2063SBarry Smith 
20185615d1e5SSatish Balay #undef __FUNC__
20195615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
202019bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
202117ab2063SBarry Smith {
2022416022c9SBarry Smith   Mat_SeqAIJ   *a;
2023416022c9SBarry Smith   Mat          B;
202417699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
2025bcd2baecSBarry Smith   MPI_Comm     comm;
202617ab2063SBarry Smith 
20273a40ed3dSBarry Smith   PetscFunctionBegin;
202819bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
202917699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
2030a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor");
203190ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
20320752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
2033a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file");
203417ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
203517ab2063SBarry Smith 
2036d64ed03dSBarry Smith   if (nz < 0) {
2037a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ");
2038d64ed03dSBarry Smith   }
2039d64ed03dSBarry Smith 
204017ab2063SBarry Smith   /* read in row lengths */
20410452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
20420752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
204317ab2063SBarry Smith 
204417ab2063SBarry Smith   /* create our matrix */
2045416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
2046416022c9SBarry Smith   B = *A;
2047416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
2048416022c9SBarry Smith   shift = a->indexshift;
204917ab2063SBarry Smith 
205017ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
20510752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr);
205217ab2063SBarry Smith   if (shift) {
205317ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
2054416022c9SBarry Smith       a->j[i] += 1;
205517ab2063SBarry Smith     }
205617ab2063SBarry Smith   }
205717ab2063SBarry Smith 
205817ab2063SBarry Smith   /* read in nonzero values */
20590752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr);
206017ab2063SBarry Smith 
206117ab2063SBarry Smith   /* set matrix "i" values */
2062416022c9SBarry Smith   a->i[0] = -shift;
206317ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
2064416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2065416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
206617ab2063SBarry Smith   }
20670452661fSBarry Smith   PetscFree(rowlengths);
206817ab2063SBarry Smith 
20696d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20706d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20713a40ed3dSBarry Smith   PetscFunctionReturn(0);
207217ab2063SBarry Smith }
207317ab2063SBarry Smith 
20745615d1e5SSatish Balay #undef __FUNC__
20755615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
20768f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
20777264ac53SSatish Balay {
20787264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
20797264ac53SSatish Balay 
20803a40ed3dSBarry Smith   PetscFunctionBegin;
2081a8c6a408SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
20827264ac53SSatish Balay 
20837264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
20847264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
2085bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
20863a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2087bcd2baecSBarry Smith   }
20887264ac53SSatish Balay 
20897264ac53SSatish Balay   /* if the a->i are the same */
20908108c231SLois Curfman McInnes   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
20913a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
20927264ac53SSatish Balay   }
20937264ac53SSatish Balay 
20947264ac53SSatish Balay   /* if a->j are the same */
2095bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
20963a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2097bcd2baecSBarry Smith   }
2098bcd2baecSBarry Smith 
2099bcd2baecSBarry Smith   /* if a->a are the same */
210019bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
21013a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
21027264ac53SSatish Balay   }
210377c4ece6SBarry Smith   *flg = PETSC_TRUE;
21043a40ed3dSBarry Smith   PetscFunctionReturn(0);
21057264ac53SSatish Balay 
21067264ac53SSatish Balay }
2107