xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 35b0346b5b6dc161b4996b5ca569fcab01f523ee)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*35b0346bSBarry Smith static char vcid[] = "$Id: aij.c,v 1.280 1998/08/03 14:58:58 balay 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"
16eeb667a2SSatish Balay #include "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++ ) {
915*35b0346bSBarry Smith     diag[i] = a->i[i+1];
916416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
917416022c9SBarry Smith       if (a->j[j]+shift == i) {
91817ab2063SBarry Smith         diag[i] = j - shift;
91917ab2063SBarry Smith         break;
92017ab2063SBarry Smith       }
92117ab2063SBarry Smith     }
92217ab2063SBarry Smith   }
923416022c9SBarry Smith   a->diag = diag;
9243a40ed3dSBarry Smith   PetscFunctionReturn(0);
92517ab2063SBarry Smith }
92617ab2063SBarry Smith 
9275615d1e5SSatish Balay #undef __FUNC__
9285615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
92944cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
93017ab2063SBarry Smith                            double fshift,int its,Vec xx)
93117ab2063SBarry Smith {
932416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
933416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
934d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
93517ab2063SBarry Smith 
9363a40ed3dSBarry Smith   PetscFunctionBegin;
937e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
938e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
939d00d2cf4SBarry Smith   if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);}
940416022c9SBarry Smith   diag = a->diag;
941416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
94217ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
94317ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
94417ab2063SBarry Smith     bs = b + shift;
94517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
946416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
947416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
948488ecbafSBarry Smith 	PLogFlops(2*n-1);
949416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
950416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
95117ab2063SBarry Smith         sum  = b[i]*d/omega;
95217ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
95317ab2063SBarry Smith         x[i] = sum;
95417ab2063SBarry Smith     }
9553a40ed3dSBarry Smith     PetscFunctionReturn(0);
95617ab2063SBarry Smith   }
95717ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
958a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done");
9593a40ed3dSBarry Smith   } else if (flag & SOR_EISENSTAT) {
96017ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
96117ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
96217ab2063SBarry Smith 
96317ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
96417ab2063SBarry Smith 
96517ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
96617ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
96717ab2063SBarry Smith     is the relaxation factor.
96817ab2063SBarry Smith     */
9690452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
97017ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
97117ab2063SBarry Smith 
97217ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
97317ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
974416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
975416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
976488ecbafSBarry Smith       PLogFlops(2*n-1);
977416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
978416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
97917ab2063SBarry Smith       sum  = b[i];
98017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
98117ab2063SBarry Smith       x[i] = omega*(sum/d);
98217ab2063SBarry Smith     }
98317ab2063SBarry Smith 
98417ab2063SBarry Smith     /*  t = b - (2*E - D)x */
985416022c9SBarry Smith     v = a->a;
98617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
987488ecbafSBarry Smith     PLogFlops(2*m);
988488ecbafSBarry Smith 
98917ab2063SBarry Smith 
99017ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
991416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
992416022c9SBarry Smith     diag = a->diag;
99317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
994416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
995416022c9SBarry Smith       n    = diag[i] - a->i[i];
996488ecbafSBarry Smith       PLogFlops(2*n-1);
997416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
998416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
99917ab2063SBarry Smith       sum  = t[i];
100017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
100117ab2063SBarry Smith       t[i] = omega*(sum/d);
100217ab2063SBarry Smith     }
100317ab2063SBarry Smith 
100417ab2063SBarry Smith     /*  x = x + t */
100517ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
10060452661fSBarry Smith     PetscFree(t);
10073a40ed3dSBarry Smith     PetscFunctionReturn(0);
100817ab2063SBarry Smith   }
100917ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
101017ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
101117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1012416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1013416022c9SBarry Smith         n    = diag[i] - a->i[i];
1014488ecbafSBarry Smith 	PLogFlops(2*n-1);
1015416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1016416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
101717ab2063SBarry Smith         sum  = b[i];
101817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
101917ab2063SBarry Smith         x[i] = omega*(sum/d);
102017ab2063SBarry Smith       }
102117ab2063SBarry Smith       xb = x;
10223a40ed3dSBarry Smith     } else xb = b;
102317ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
102417ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
102517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1026416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
102717ab2063SBarry Smith       }
1028488ecbafSBarry Smith       PLogFlops(m);
102917ab2063SBarry Smith     }
103017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
103117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1032416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1033416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
1034488ecbafSBarry Smith 	PLogFlops(2*n-1);
1035416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1036416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
103717ab2063SBarry Smith         sum  = xb[i];
103817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
103917ab2063SBarry Smith         x[i] = omega*(sum/d);
104017ab2063SBarry Smith       }
104117ab2063SBarry Smith     }
104217ab2063SBarry Smith     its--;
104317ab2063SBarry Smith   }
104417ab2063SBarry Smith   while (its--) {
104517ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
104617ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1047416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1048416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1049488ecbafSBarry Smith 	PLogFlops(2*n-1);
1050416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1051416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
105217ab2063SBarry Smith         sum  = b[i];
105317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10547e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
105517ab2063SBarry Smith       }
105617ab2063SBarry Smith     }
105717ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
105817ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1059416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1060416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1061488ecbafSBarry Smith 	PLogFlops(2*n-1);
1062416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1063416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
106417ab2063SBarry Smith         sum  = b[i];
106517ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10667e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
106717ab2063SBarry Smith       }
106817ab2063SBarry Smith     }
106917ab2063SBarry Smith   }
10703a40ed3dSBarry Smith   PetscFunctionReturn(0);
107117ab2063SBarry Smith }
107217ab2063SBarry Smith 
10735615d1e5SSatish Balay #undef __FUNC__
10745615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10758f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
107617ab2063SBarry Smith {
1077416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
10784e220ebcSLois Curfman McInnes 
10793a40ed3dSBarry Smith   PetscFunctionBegin;
10804e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
10814e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
10824e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
10834e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
10844e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10854e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
10864e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
10874e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
10884e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
10894e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
10904e220ebcSLois Curfman McInnes   info->memory         = A->mem;
10914e220ebcSLois Curfman McInnes   if (A->factor) {
10924e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10934e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
10944e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
10954e220ebcSLois Curfman McInnes   } else {
10964e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
10974e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
10984e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
10994e220ebcSLois Curfman McInnes   }
11003a40ed3dSBarry Smith   PetscFunctionReturn(0);
110117ab2063SBarry Smith }
110217ab2063SBarry Smith 
110317ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
110417ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
110517ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
110617ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
110717ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
110817ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
110917ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
111017ab2063SBarry Smith 
11115615d1e5SSatish Balay #undef __FUNC__
11125615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
11138f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
111417ab2063SBarry Smith {
1115416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1116416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
111717ab2063SBarry Smith 
11183a40ed3dSBarry Smith   PetscFunctionBegin;
111977c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
112017ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
112117ab2063SBarry Smith   if (diag) {
112217ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1123a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1124416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1125416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1126416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1127416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
11283a40ed3dSBarry Smith       } else {
1129d64ed03dSBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr);
113017ab2063SBarry Smith       }
113117ab2063SBarry Smith     }
11323a40ed3dSBarry Smith   } else {
113317ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1134a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1135416022c9SBarry Smith       a->ilen[rows[i]] = 0;
113617ab2063SBarry Smith     }
113717ab2063SBarry Smith   }
113817ab2063SBarry Smith   ISRestoreIndices(is,&rows);
113943a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11403a40ed3dSBarry Smith   PetscFunctionReturn(0);
114117ab2063SBarry Smith }
114217ab2063SBarry Smith 
11435615d1e5SSatish Balay #undef __FUNC__
11445615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11458f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
114617ab2063SBarry Smith {
1147416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11483a40ed3dSBarry Smith 
11493a40ed3dSBarry Smith   PetscFunctionBegin;
11500752156aSBarry Smith   if (m) *m = a->m;
11510752156aSBarry Smith   if (n) *n = a->n;
11523a40ed3dSBarry Smith   PetscFunctionReturn(0);
115317ab2063SBarry Smith }
115417ab2063SBarry Smith 
11555615d1e5SSatish Balay #undef __FUNC__
11565615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11578f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
115817ab2063SBarry Smith {
1159416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11603a40ed3dSBarry Smith 
11613a40ed3dSBarry Smith   PetscFunctionBegin;
1162416022c9SBarry Smith   *m = 0; *n = a->m;
11633a40ed3dSBarry Smith   PetscFunctionReturn(0);
116417ab2063SBarry Smith }
1165026e39d0SSatish Balay 
11665615d1e5SSatish Balay #undef __FUNC__
11675615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11684e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
116917ab2063SBarry Smith {
1170416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1171c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
117217ab2063SBarry Smith 
11733a40ed3dSBarry Smith   PetscFunctionBegin;
1174a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
117517ab2063SBarry Smith 
1176416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1177416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
117817ab2063SBarry Smith   if (idx) {
1179416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
11804e093b46SBarry Smith     if (*nz && shift) {
11810452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
118217ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
11834e093b46SBarry Smith     } else if (*nz) {
11844e093b46SBarry Smith       *idx = itmp;
118517ab2063SBarry Smith     }
118617ab2063SBarry Smith     else *idx = 0;
118717ab2063SBarry Smith   }
11883a40ed3dSBarry Smith   PetscFunctionReturn(0);
118917ab2063SBarry Smith }
119017ab2063SBarry Smith 
11915615d1e5SSatish Balay #undef __FUNC__
11925615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
11934e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
119417ab2063SBarry Smith {
11954e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11963a40ed3dSBarry Smith 
11973a40ed3dSBarry Smith   PetscFunctionBegin;
11984e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
11993a40ed3dSBarry Smith   PetscFunctionReturn(0);
120017ab2063SBarry Smith }
120117ab2063SBarry Smith 
12025615d1e5SSatish Balay #undef __FUNC__
12035615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
12048f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
120517ab2063SBarry Smith {
1206416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1207416022c9SBarry Smith   Scalar     *v = a->a;
120817ab2063SBarry Smith   double     sum = 0.0;
1209416022c9SBarry Smith   int        i, j,shift = a->indexshift;
121017ab2063SBarry Smith 
12113a40ed3dSBarry Smith   PetscFunctionBegin;
121217ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1213416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
12143a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
1215e20fef11SSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
121617ab2063SBarry Smith #else
121717ab2063SBarry Smith       sum += (*v)*(*v); v++;
121817ab2063SBarry Smith #endif
121917ab2063SBarry Smith     }
122017ab2063SBarry Smith     *norm = sqrt(sum);
12213a40ed3dSBarry Smith   } else if (type == NORM_1) {
122217ab2063SBarry Smith     double *tmp;
1223416022c9SBarry Smith     int    *jj = a->j;
122466963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1225cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
122617ab2063SBarry Smith     *norm = 0.0;
1227416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1228a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
122917ab2063SBarry Smith     }
1230416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
123117ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
123217ab2063SBarry Smith     }
12330452661fSBarry Smith     PetscFree(tmp);
12343a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
123517ab2063SBarry Smith     *norm = 0.0;
1236416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1237416022c9SBarry Smith       v = a->a + a->i[j] + shift;
123817ab2063SBarry Smith       sum = 0.0;
1239416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1240cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
124117ab2063SBarry Smith       }
124217ab2063SBarry Smith       if (sum > *norm) *norm = sum;
124317ab2063SBarry Smith     }
12443a40ed3dSBarry Smith   } else {
1245a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
124617ab2063SBarry Smith   }
12473a40ed3dSBarry Smith   PetscFunctionReturn(0);
124817ab2063SBarry Smith }
124917ab2063SBarry Smith 
12505615d1e5SSatish Balay #undef __FUNC__
12515615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12528f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
125317ab2063SBarry Smith {
1254416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1255416022c9SBarry Smith   Mat        C;
1256416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1257416022c9SBarry Smith   int        shift = a->indexshift;
1258d5d45c9bSBarry Smith   Scalar     *array = a->a;
125917ab2063SBarry Smith 
12603a40ed3dSBarry Smith   PetscFunctionBegin;
1261a8c6a408SBarry Smith   if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place");
12620452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1263cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
126417ab2063SBarry Smith   if (shift) {
126517ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
126617ab2063SBarry Smith   }
126717ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1268416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12690452661fSBarry Smith   PetscFree(col);
127017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
127117ab2063SBarry Smith     len = ai[i+1]-ai[i];
1272416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
127317ab2063SBarry Smith     array += len; aj += len;
127417ab2063SBarry Smith   }
127517ab2063SBarry Smith   if (shift) {
127617ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
127717ab2063SBarry Smith   }
127817ab2063SBarry Smith 
12796d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12806d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
128117ab2063SBarry Smith 
12823638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1283416022c9SBarry Smith     *B = C;
128417ab2063SBarry Smith   } else {
1285f830108cSBarry Smith     PetscOps *Abops;
12860a6ffc59SBarry Smith     MatOps   Aops;
1287f830108cSBarry Smith 
1288416022c9SBarry Smith     /* This isn't really an in-place transpose */
12890452661fSBarry Smith     PetscFree(a->a);
12900452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
12910452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
12920452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
12930452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
12940452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
12951073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
12960452661fSBarry Smith     PetscFree(a);
1297f830108cSBarry Smith 
1298488ecbafSBarry Smith 
1299488ecbafSBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
1300488ecbafSBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
1301488ecbafSBarry Smith 
1302f830108cSBarry Smith     /*
1303f830108cSBarry Smith       This is horrible, horrible code. We need to keep the
13048f0f457cSSatish Balay       the bops and ops Structures, copy everything from C
13058f0f457cSSatish Balay       including the function pointers..
1306f830108cSBarry Smith     */
13078f0f457cSSatish Balay     PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps));
13088f0f457cSSatish Balay     PetscMemcpy(A->bops,C->bops,sizeof(PetscOps));
1309f830108cSBarry Smith     Abops    = A->bops;
1310f830108cSBarry Smith     Aops     = A->ops;
1311f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
1312f830108cSBarry Smith     A->bops  = Abops;
1313f830108cSBarry Smith     A->ops   = Aops;
131427a8da17SBarry Smith     A->qlist = 0;
1315f830108cSBarry Smith 
13160452661fSBarry Smith     PetscHeaderDestroy(C);
131717ab2063SBarry Smith   }
13183a40ed3dSBarry Smith   PetscFunctionReturn(0);
131917ab2063SBarry Smith }
132017ab2063SBarry Smith 
13215615d1e5SSatish Balay #undef __FUNC__
13225615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
13238f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
132417ab2063SBarry Smith {
1325416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
132617ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1327e1311b90SBarry Smith   int        ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
132817ab2063SBarry Smith 
13293a40ed3dSBarry Smith   PetscFunctionBegin;
133017ab2063SBarry Smith   if (ll) {
13313ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
13323ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
1333e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr);
1334a8c6a408SBarry Smith     if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
1335e1311b90SBarry Smith     ierr = VecGetArray(ll,&l); CHKERRQ(ierr);
1336416022c9SBarry Smith     v = a->a;
133717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
133817ab2063SBarry Smith       x = l[i];
1339416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
134017ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
134117ab2063SBarry Smith     }
1342e1311b90SBarry Smith     ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr);
134344cd7ae7SLois Curfman McInnes     PLogFlops(nz);
134417ab2063SBarry Smith   }
134517ab2063SBarry Smith   if (rr) {
1346e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr);
1347a8c6a408SBarry Smith     if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
1348e1311b90SBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1349416022c9SBarry Smith     v = a->a; jj = a->j;
135017ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
135117ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
135217ab2063SBarry Smith     }
1353e1311b90SBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
135444cd7ae7SLois Curfman McInnes     PLogFlops(nz);
135517ab2063SBarry Smith   }
13563a40ed3dSBarry Smith   PetscFunctionReturn(0);
135717ab2063SBarry Smith }
135817ab2063SBarry Smith 
13595615d1e5SSatish Balay #undef __FUNC__
13605615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
13616a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B)
136217ab2063SBarry Smith {
1363db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1364d5db1b72SBarry Smith   int          *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
136599141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1366a2744918SBarry Smith   register int sum,lensi;
136799141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
136899141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
136999141d43SSatish Balay   Scalar       *a_new,*mat_a;
1370416022c9SBarry Smith   Mat          C;
1371fee21e36SBarry Smith   PetscTruth   stride;
137217ab2063SBarry Smith 
13733a40ed3dSBarry Smith   PetscFunctionBegin;
1374d64ed03dSBarry Smith   ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr);
1375a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
1376d64ed03dSBarry Smith   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
1377a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
137899141d43SSatish Balay 
137917ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
138017ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
138117ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
138217ab2063SBarry Smith 
1383fee21e36SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr);
1384fee21e36SBarry Smith   ierr = ISStride(iscol,&stride); CHKERRQ(ierr);
1385fee21e36SBarry Smith   if (stride && step == 1) {
138602834360SBarry Smith     /* special case of contiguous rows */
138757faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
138802834360SBarry Smith     starts = lens + ncols;
138902834360SBarry Smith     /* loop over new rows determining lens and starting points */
139002834360SBarry Smith     for (i=0; i<nrows; i++) {
1391a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1392a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
139302834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1394d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
139502834360SBarry Smith           starts[i] = k;
139602834360SBarry Smith           break;
139702834360SBarry Smith 	}
139802834360SBarry Smith       }
1399a2744918SBarry Smith       sum = 0;
140002834360SBarry Smith       while (k < kend) {
1401d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1402a2744918SBarry Smith         sum++;
140302834360SBarry Smith       }
1404a2744918SBarry Smith       lens[i] = sum;
140502834360SBarry Smith     }
140602834360SBarry Smith     /* create submatrix */
1407cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
140808480c60SBarry Smith       int n_cols,n_rows;
140908480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1410a8c6a408SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1411d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
141208480c60SBarry Smith       C = *B;
14133a40ed3dSBarry Smith     } else {
141402834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
141508480c60SBarry Smith     }
1416db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1417db02288aSLois Curfman McInnes 
141802834360SBarry Smith     /* loop over rows inserting into submatrix */
1419db02288aSLois Curfman McInnes     a_new    = c->a;
1420db02288aSLois Curfman McInnes     j_new    = c->j;
1421db02288aSLois Curfman McInnes     i_new    = c->i;
1422db02288aSLois Curfman McInnes     i_new[0] = -shift;
142302834360SBarry Smith     for (i=0; i<nrows; i++) {
1424a2744918SBarry Smith       ii    = starts[i];
1425a2744918SBarry Smith       lensi = lens[i];
1426a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1427a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
142802834360SBarry Smith       }
1429a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1430a2744918SBarry Smith       a_new      += lensi;
1431a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1432a2744918SBarry Smith       c->ilen[i]  = lensi;
143302834360SBarry Smith     }
14340452661fSBarry Smith     PetscFree(lens);
14353a40ed3dSBarry Smith   } else {
143602834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
14370452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
143802834360SBarry Smith     ssmap = smap + shift;
143999141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1440cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
144117ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
144202834360SBarry Smith     /* determine lens of each row */
144302834360SBarry Smith     for (i=0; i<nrows; i++) {
1444d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
144502834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
144602834360SBarry Smith       lens[i] = 0;
144702834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1448d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
144902834360SBarry Smith           lens[i]++;
145002834360SBarry Smith         }
145102834360SBarry Smith       }
145202834360SBarry Smith     }
145317ab2063SBarry Smith     /* Create and fill new matrix */
1454a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
145599141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
1456a8c6a408SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
145799141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1458a8c6a408SBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
145999141d43SSatish Balay       }
146099141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
146108480c60SBarry Smith       C = *B;
14623a40ed3dSBarry Smith     } else {
146302834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
146408480c60SBarry Smith     }
146599141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
146617ab2063SBarry Smith     for (i=0; i<nrows; i++) {
146799141d43SSatish Balay       row    = irow[i];
146899141d43SSatish Balay       kstart = ai[row]+shift;
146999141d43SSatish Balay       kend   = kstart + a->ilen[row];
147099141d43SSatish Balay       mat_i  = c->i[i]+shift;
147199141d43SSatish Balay       mat_j  = c->j + mat_i;
147299141d43SSatish Balay       mat_a  = c->a + mat_i;
147399141d43SSatish Balay       mat_ilen = c->ilen + i;
147417ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
147599141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
147699141d43SSatish Balay           *mat_j++ = tcol - (!shift);
147799141d43SSatish Balay           *mat_a++ = a->a[k];
147899141d43SSatish Balay           (*mat_ilen)++;
147999141d43SSatish Balay 
148017ab2063SBarry Smith         }
148117ab2063SBarry Smith       }
148217ab2063SBarry Smith     }
148302834360SBarry Smith     /* Free work space */
148402834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
148599141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
148602834360SBarry Smith   }
14876d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
14886d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
148917ab2063SBarry Smith 
149017ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1491416022c9SBarry Smith   *B = C;
14923a40ed3dSBarry Smith   PetscFunctionReturn(0);
149317ab2063SBarry Smith }
149417ab2063SBarry Smith 
1495a871dcd8SBarry Smith /*
149663b91edcSBarry Smith      note: This can only work for identity for row and col. It would
149763b91edcSBarry Smith    be good to check this and otherwise generate an error.
1498a871dcd8SBarry Smith */
14995615d1e5SSatish Balay #undef __FUNC__
15005615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
15018f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1502a871dcd8SBarry Smith {
150363b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
150408480c60SBarry Smith   int        ierr;
150563b91edcSBarry Smith   Mat        outA;
150663b91edcSBarry Smith 
15073a40ed3dSBarry Smith   PetscFunctionBegin;
1508a8c6a408SBarry Smith   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
1509a871dcd8SBarry Smith 
151063b91edcSBarry Smith   outA          = inA;
151163b91edcSBarry Smith   inA->factor   = FACTOR_LU;
151263b91edcSBarry Smith   a->row        = row;
151363b91edcSBarry Smith   a->col        = col;
151463b91edcSBarry Smith 
1515f0ec6fceSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1516f0ec6fceSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
15171e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1518f0ec6fceSSatish Balay 
151994a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
15200452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work);
152194a9d846SBarry Smith   }
152263b91edcSBarry Smith 
152308480c60SBarry Smith   if (!a->diag) {
152408480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
152563b91edcSBarry Smith   }
152663b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
15273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1528a871dcd8SBarry Smith }
1529a871dcd8SBarry Smith 
153074cdf0dfSBarry Smith #include "pinclude/blaslapack.h"
15315615d1e5SSatish Balay #undef __FUNC__
15325615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
15338f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1534f0b747eeSBarry Smith {
1535f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1536f0b747eeSBarry Smith   int        one = 1;
15373a40ed3dSBarry Smith 
15383a40ed3dSBarry Smith   PetscFunctionBegin;
1539f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1540f0b747eeSBarry Smith   PLogFlops(a->nz);
15413a40ed3dSBarry Smith   PetscFunctionReturn(0);
1542f0b747eeSBarry Smith }
1543f0b747eeSBarry Smith 
15445615d1e5SSatish Balay #undef __FUNC__
15455615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
15466a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1547cddf8d76SBarry Smith {
1548cddf8d76SBarry Smith   int ierr,i;
1549cddf8d76SBarry Smith 
15503a40ed3dSBarry Smith   PetscFunctionBegin;
1551cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
15520452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1553cddf8d76SBarry Smith   }
1554cddf8d76SBarry Smith 
1555cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
15566a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1557cddf8d76SBarry Smith   }
15583a40ed3dSBarry Smith   PetscFunctionReturn(0);
1559cddf8d76SBarry Smith }
1560cddf8d76SBarry Smith 
15615615d1e5SSatish Balay #undef __FUNC__
15625615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
15638f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
15645a838052SSatish Balay {
1565f830108cSBarry Smith   PetscFunctionBegin;
15665a838052SSatish Balay   *bs = 1;
15673a40ed3dSBarry Smith   PetscFunctionReturn(0);
15685a838052SSatish Balay }
15695a838052SSatish Balay 
15705615d1e5SSatish Balay #undef __FUNC__
15715615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15728f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15734dcbc457SBarry Smith {
1574e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
157506763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15768a047759SSatish Balay   int        start, end, *ai, *aj;
1577bbd702dbSSatish Balay   BT         table;
1578bbd702dbSSatish Balay 
15793a40ed3dSBarry Smith   PetscFunctionBegin;
15808a047759SSatish Balay   shift = a->indexshift;
1581e4d965acSSatish Balay   m     = a->m;
1582e4d965acSSatish Balay   ai    = a->i;
15838a047759SSatish Balay   aj    = a->j+shift;
15848a047759SSatish Balay 
1585a8c6a408SBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used");
158606763907SSatish Balay 
158706763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1588bbd702dbSSatish Balay   ierr  = BTCreate(m,table); CHKERRQ(ierr);
158906763907SSatish Balay 
1590e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1591b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1592e4d965acSSatish Balay     isz  = 0;
1593bbd702dbSSatish Balay     BTMemzero(m,table);
1594e4d965acSSatish Balay 
1595e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
15964dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
159777c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1598e4d965acSSatish Balay 
1599dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1600e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
1601bbd702dbSSatish Balay       if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];}
16024dcbc457SBarry Smith     }
160306763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
160406763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1605e4d965acSSatish Balay 
160604a348a9SBarry Smith     k = 0;
160704a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
160804a348a9SBarry Smith       n = isz;
160906763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1610e4d965acSSatish Balay         row   = nidx[k];
1611e4d965acSSatish Balay         start = ai[row];
1612e4d965acSSatish Balay         end   = ai[row+1];
161304a348a9SBarry Smith         for ( l = start; l<end ; l++){
16148a047759SSatish Balay           val = aj[l] + shift;
16152a8f2162SSatish Balay           if (!BTLookupSet(table,val)) {nidx[isz++] = val;}
1616e4d965acSSatish Balay         }
1617e4d965acSSatish Balay       }
1618e4d965acSSatish Balay     }
1619029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1620e4d965acSSatish Balay   }
1621bbd702dbSSatish Balay   BTDestroy(table);
162206763907SSatish Balay   PetscFree(nidx);
16233a40ed3dSBarry Smith   PetscFunctionReturn(0);
16244dcbc457SBarry Smith }
162517ab2063SBarry Smith 
16260513a670SBarry Smith /* -------------------------------------------------------------- */
16270513a670SBarry Smith #undef __FUNC__
16280513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
16290513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
16300513a670SBarry Smith {
16310513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
16320513a670SBarry Smith   Scalar     *vwork;
16330513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
16340513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
163556cd22aeSBarry Smith   IS         icolp,irowp;
16360513a670SBarry Smith 
16373a40ed3dSBarry Smith   PetscFunctionBegin;
163856cd22aeSBarry Smith   ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr);
163956cd22aeSBarry Smith   ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr);
164056cd22aeSBarry Smith   ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr);
164156cd22aeSBarry Smith   ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr);
16420513a670SBarry Smith 
16430513a670SBarry Smith   /* determine lengths of permuted rows */
16440513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
16450513a670SBarry Smith   for (i=0; i<m; i++ ) {
16460513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
16470513a670SBarry Smith   }
16480513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
16490513a670SBarry Smith   PetscFree(lens);
16500513a670SBarry Smith 
16510513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
16520513a670SBarry Smith   for (i=0; i<m; i++) {
16530513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16540513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
16550513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
16560513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16570513a670SBarry Smith   }
16580513a670SBarry Smith   PetscFree(cnew);
16590513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
16600513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
166156cd22aeSBarry Smith   ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr);
166256cd22aeSBarry Smith   ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr);
166356cd22aeSBarry Smith   ierr = ISDestroy(irowp); CHKERRQ(ierr);
166456cd22aeSBarry Smith   ierr = ISDestroy(icolp); CHKERRQ(ierr);
16653a40ed3dSBarry Smith   PetscFunctionReturn(0);
16660513a670SBarry Smith }
16670513a670SBarry Smith 
16685615d1e5SSatish Balay #undef __FUNC__
16695615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1670682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1671682d7d0cSBarry Smith {
1672682d7d0cSBarry Smith   static int called = 0;
1673682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1674682d7d0cSBarry Smith 
16753a40ed3dSBarry Smith   PetscFunctionBegin;
16763a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
167776be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
167876be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
167976be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
168076be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodes\n");
168176be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1682682d7d0cSBarry Smith #if defined(HAVE_ESSL)
168376be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1684682d7d0cSBarry Smith #endif
16853a40ed3dSBarry Smith   PetscFunctionReturn(0);
1686682d7d0cSBarry Smith }
16878f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1688a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1689a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1690a93ec695SBarry Smith 
1691682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
16920a6ffc59SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
169317ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1694416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1695416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
169617ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
169717ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
169817ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
169917ab2063SBarry Smith        MatRelax_SeqAIJ,
170017ab2063SBarry Smith        MatTranspose_SeqAIJ,
17017264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1702f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
170317ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
170417ab2063SBarry Smith        MatCompress_SeqAIJ,
170517ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
170617ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
170717ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
170817ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
170994a9d846SBarry Smith        0,0,
17103d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1711cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
17127eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1713682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1714f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
17155a838052SSatish Balay        MatScale_SeqAIJ,0,0,
17166945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
17176945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
17183b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
17193b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
17203b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1721a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1722a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
17230513a670SBarry Smith        MatColoringPatch_SeqAIJ,
17240513a670SBarry Smith        0,
1725cda55fadSBarry Smith        MatPermute_SeqAIJ,
1726cda55fadSBarry Smith        0,
1727cda55fadSBarry Smith        0,
1728cda55fadSBarry Smith        0,
1729cda55fadSBarry Smith        0,
1730cda55fadSBarry Smith        MatGetMaps_Petsc};
173117ab2063SBarry Smith 
173217ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
173317ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
173417ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
173517ab2063SBarry Smith 
17365615d1e5SSatish Balay #undef __FUNC__
1737bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ"
1738bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
1739bef8e0ddSBarry Smith {
1740bef8e0ddSBarry Smith   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
1741bef8e0ddSBarry Smith   int        i,nz,n;
1742bef8e0ddSBarry Smith 
1743bef8e0ddSBarry Smith   PetscFunctionBegin;
1744bef8e0ddSBarry Smith   if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing");
1745bef8e0ddSBarry Smith 
1746bef8e0ddSBarry Smith   nz = aij->maxnz;
1747bef8e0ddSBarry Smith   n  = aij->n;
1748bef8e0ddSBarry Smith   for (i=0; i<nz; i++) {
1749bef8e0ddSBarry Smith     aij->j[i] = indices[i];
1750bef8e0ddSBarry Smith   }
1751bef8e0ddSBarry Smith   aij->nz = nz;
1752bef8e0ddSBarry Smith   for ( i=0; i<n; i++ ) {
1753bef8e0ddSBarry Smith     aij->ilen[i] = aij->imax[i];
1754bef8e0ddSBarry Smith   }
1755bef8e0ddSBarry Smith 
1756bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1757bef8e0ddSBarry Smith }
1758bef8e0ddSBarry Smith 
1759bef8e0ddSBarry Smith #undef __FUNC__
1760bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices"
1761bef8e0ddSBarry Smith /*@
1762bef8e0ddSBarry Smith     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
1763bef8e0ddSBarry Smith        in the matrix.
1764bef8e0ddSBarry Smith 
1765bef8e0ddSBarry Smith   Input Parameters:
1766bef8e0ddSBarry Smith +  mat - the SeqAIJ matrix
1767bef8e0ddSBarry Smith -  indices - the column indices
1768bef8e0ddSBarry Smith 
1769bef8e0ddSBarry Smith   Notes:
1770bef8e0ddSBarry Smith     This can be called if you have precomputed the nonzero structure of the
1771bef8e0ddSBarry Smith   matrix and want to provide it to the matrix object to improve the performance
1772bef8e0ddSBarry Smith   of the MatSetValues() operation.
1773bef8e0ddSBarry Smith 
1774bef8e0ddSBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
1775bef8e0ddSBarry Smith   MatCreateSeqAIJ().
1776bef8e0ddSBarry Smith 
1777bef8e0ddSBarry Smith     MUST be called before any calls to MatSetValues();
1778bef8e0ddSBarry Smith 
1779bef8e0ddSBarry Smith @*/
1780bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
1781bef8e0ddSBarry Smith {
1782bef8e0ddSBarry Smith   int ierr,(*f)(Mat,int *);
1783bef8e0ddSBarry Smith 
1784bef8e0ddSBarry Smith   PetscFunctionBegin;
1785bef8e0ddSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1786bef8e0ddSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f);
1787bef8e0ddSBarry Smith          CHKERRQ(ierr);
1788bef8e0ddSBarry Smith   if (f) {
1789bef8e0ddSBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
1790bef8e0ddSBarry Smith   } else {
1791bef8e0ddSBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
1792bef8e0ddSBarry Smith   }
1793bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1794bef8e0ddSBarry Smith }
1795bef8e0ddSBarry Smith 
1796bef8e0ddSBarry Smith #undef __FUNC__
17975615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
179817ab2063SBarry Smith /*@C
1799682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
18000d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
18016e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
18022bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
18032bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
180417ab2063SBarry Smith 
1805db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1806db81eaa0SLois Curfman McInnes 
180717ab2063SBarry Smith    Input Parameters:
1808db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
180917ab2063SBarry Smith .  m - number of rows
181017ab2063SBarry Smith .  n - number of columns
181117ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
1812db81eaa0SLois Curfman McInnes -  nzz - array containing the number of nonzeros in the various rows
18132bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
181417ab2063SBarry Smith 
181517ab2063SBarry Smith    Output Parameter:
1816416022c9SBarry Smith .  A - the matrix
181717ab2063SBarry Smith 
1818b259b22eSLois Curfman McInnes    Notes:
181917ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
182017ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
18210002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
182244cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
182317ab2063SBarry Smith 
182417ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1825a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
18263d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
18276da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
182817ab2063SBarry Smith 
1829682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
18304fca80b9SLois Curfman McInnes    improve numerical efficiency of matrix-vector products and solves. We
1831682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
18326c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
18336c7ebb05SLois Curfman McInnes 
18346c7ebb05SLois Curfman McInnes    Options Database Keys:
1835db81eaa0SLois Curfman McInnes +  -mat_aij_no_inode  - Do not use inodes
1836db81eaa0SLois Curfman McInnes .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
1837db81eaa0SLois Curfman McInnes -  -mat_aij_oneindex - Internally use indexing starting at 1
1838db81eaa0SLois Curfman McInnes         rather than 0.  Note that when calling MatSetValues(),
1839db81eaa0SLois Curfman McInnes         the user still MUST index entries starting at 0!
184017ab2063SBarry Smith 
1841bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices()
184217ab2063SBarry Smith @*/
1843416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
184417ab2063SBarry Smith {
1845416022c9SBarry Smith   Mat        B;
1846416022c9SBarry Smith   Mat_SeqAIJ *b;
18476945ee14SBarry Smith   int        i, len, ierr, flg,size;
18486945ee14SBarry Smith 
18493a40ed3dSBarry Smith   PetscFunctionBegin;
18506945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1851a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1");
1852d5d45c9bSBarry Smith 
1853416022c9SBarry Smith   *A                  = 0;
1854f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView);
1855416022c9SBarry Smith   PLogObjectCreate(B);
18560452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
185744cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
18580a6ffc59SBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1859e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_SeqAIJ;
1860e1311b90SBarry Smith   B->ops->view             = MatView_SeqAIJ;
1861416022c9SBarry Smith   B->factor           = 0;
1862416022c9SBarry Smith   B->lupivotthreshold = 1.0;
186390f02eecSBarry Smith   B->mapping          = 0;
1864e1311b90SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg);CHKERRQ(ierr);
18657a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
1866e1311b90SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*)&b->ilu_preserve_row_sums);CHKERRQ(ierr);
1867416022c9SBarry Smith   b->row              = 0;
1868416022c9SBarry Smith   b->col              = 0;
186982bf6240SBarry Smith   b->icol             = 0;
1870416022c9SBarry Smith   b->indexshift       = 0;
1871b810aeb4SBarry Smith   b->reallocs         = 0;
187269957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
187369957df2SSatish Balay   if (flg) b->indexshift = -1;
187417ab2063SBarry Smith 
187544cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
187644cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
1877a5ae1ecdSBarry Smith 
18784d197716SBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
18794d197716SBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1880a5ae1ecdSBarry Smith 
18810452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1882b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
18837b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
18847b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1885416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
188617ab2063SBarry Smith     nz = nz*m;
18873a40ed3dSBarry Smith   } else {
188817ab2063SBarry Smith     nz = 0;
1889416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
189017ab2063SBarry Smith   }
189117ab2063SBarry Smith 
189217ab2063SBarry Smith   /* allocate the matrix space */
1893416022c9SBarry Smith   len             = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
18940452661fSBarry Smith   b->a            = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1895416022c9SBarry Smith   b->j            = (int *) (b->a + nz);
1896cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1897416022c9SBarry Smith   b->i            = b->j + nz;
1898416022c9SBarry Smith   b->singlemalloc = 1;
189917ab2063SBarry Smith 
1900416022c9SBarry Smith   b->i[0] = -b->indexshift;
190117ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1902416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
190317ab2063SBarry Smith   }
190417ab2063SBarry Smith 
1905416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
19060452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1907f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1908416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
190917ab2063SBarry Smith 
1910416022c9SBarry Smith   b->nz               = 0;
1911416022c9SBarry Smith   b->maxnz            = nz;
1912416022c9SBarry Smith   b->sorted           = 0;
1913416022c9SBarry Smith   b->roworiented      = 1;
1914416022c9SBarry Smith   b->nonew            = 0;
1915416022c9SBarry Smith   b->diag             = 0;
1916416022c9SBarry Smith   b->solve_work       = 0;
191771bd300dSLois Curfman McInnes   b->spptr            = 0;
1918754ec7b1SSatish Balay   b->inode.node_count = 0;
1919754ec7b1SSatish Balay   b->inode.size       = 0;
19206c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
19216c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
19224e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
192317ab2063SBarry Smith 
1924416022c9SBarry Smith   *A = B;
19254e220ebcSLois Curfman McInnes 
19264b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
19274b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
192869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
192969957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
19304b14c69eSBarry Smith #endif
193169957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
193269957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
193369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
193469957df2SSatish Balay   if (flg) {
1935a8c6a408SBarry Smith     if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1936416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
193717ab2063SBarry Smith   }
193869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
193969957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1940bef8e0ddSBarry Smith 
1941bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
1942bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
1943bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
19443a40ed3dSBarry Smith   PetscFunctionReturn(0);
194517ab2063SBarry Smith }
194617ab2063SBarry Smith 
19475615d1e5SSatish Balay #undef __FUNC__
19485615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
19493d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
195017ab2063SBarry Smith {
1951416022c9SBarry Smith   Mat        C;
1952416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1953bef8e0ddSBarry Smith   int        i,len, m = a->m,shift = a->indexshift,ierr;
195417ab2063SBarry Smith 
19553a40ed3dSBarry Smith   PetscFunctionBegin;
19564043dd9cSLois Curfman McInnes   *B = 0;
1957f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView);
1958416022c9SBarry Smith   PLogObjectCreate(C);
19590452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1960f830108cSBarry Smith   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1961e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqAIJ;
1962e1311b90SBarry Smith   C->ops->view       = MatView_SeqAIJ;
1963416022c9SBarry Smith   C->factor     = A->factor;
1964416022c9SBarry Smith   c->row        = 0;
1965416022c9SBarry Smith   c->col        = 0;
196682bf6240SBarry Smith   c->icol       = 0;
1967416022c9SBarry Smith   c->indexshift = shift;
1968c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
196917ab2063SBarry Smith 
197044cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
197144cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
197244cd7ae7SLois Curfman McInnes   C->M          = a->m;
197344cd7ae7SLois Curfman McInnes   C->N          = a->n;
197417ab2063SBarry Smith 
19750452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
19760452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
197717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1978416022c9SBarry Smith     c->imax[i] = a->imax[i];
1979416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
198017ab2063SBarry Smith   }
198117ab2063SBarry Smith 
198217ab2063SBarry Smith   /* allocate the matrix space */
1983416022c9SBarry Smith   c->singlemalloc = 1;
1984416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
19850452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1986416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1987416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1988416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
198917ab2063SBarry Smith   if (m > 0) {
1990416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
199108480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1992416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
199317ab2063SBarry Smith     }
199408480c60SBarry Smith   }
199517ab2063SBarry Smith 
1996f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1997416022c9SBarry Smith   c->sorted      = a->sorted;
1998416022c9SBarry Smith   c->roworiented = a->roworiented;
1999416022c9SBarry Smith   c->nonew       = a->nonew;
20007a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
200117ab2063SBarry Smith 
2002416022c9SBarry Smith   if (a->diag) {
20030452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
2004416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
200517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
2006416022c9SBarry Smith       c->diag[i] = a->diag[i];
200717ab2063SBarry Smith     }
20083a40ed3dSBarry Smith   } else c->diag          = 0;
20096c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
20106c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
2011754ec7b1SSatish Balay   if (a->inode.size){
2012daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
2013754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
2014daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
2015754ec7b1SSatish Balay   } else {
2016754ec7b1SSatish Balay     c->inode.size       = 0;
2017754ec7b1SSatish Balay     c->inode.node_count = 0;
2018754ec7b1SSatish Balay   }
2019416022c9SBarry Smith   c->nz                 = a->nz;
2020416022c9SBarry Smith   c->maxnz              = a->maxnz;
2021416022c9SBarry Smith   c->solve_work         = 0;
202276dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
2023754ec7b1SSatish Balay 
2024416022c9SBarry Smith   *B = C;
2025bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqAIJSetColumnIndices_C",
2026bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2027bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
20283a40ed3dSBarry Smith   PetscFunctionReturn(0);
202917ab2063SBarry Smith }
203017ab2063SBarry Smith 
20315615d1e5SSatish Balay #undef __FUNC__
20325615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
203319bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
203417ab2063SBarry Smith {
2035416022c9SBarry Smith   Mat_SeqAIJ   *a;
2036416022c9SBarry Smith   Mat          B;
203717699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
2038bcd2baecSBarry Smith   MPI_Comm     comm;
203917ab2063SBarry Smith 
20403a40ed3dSBarry Smith   PetscFunctionBegin;
204119bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
204217699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
2043a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor");
204490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
20450752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
2046a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file");
204717ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
204817ab2063SBarry Smith 
2049d64ed03dSBarry Smith   if (nz < 0) {
2050a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ");
2051d64ed03dSBarry Smith   }
2052d64ed03dSBarry Smith 
205317ab2063SBarry Smith   /* read in row lengths */
20540452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
20550752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
205617ab2063SBarry Smith 
205717ab2063SBarry Smith   /* create our matrix */
2058416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
2059416022c9SBarry Smith   B = *A;
2060416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
2061416022c9SBarry Smith   shift = a->indexshift;
206217ab2063SBarry Smith 
206317ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
20640752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr);
206517ab2063SBarry Smith   if (shift) {
206617ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
2067416022c9SBarry Smith       a->j[i] += 1;
206817ab2063SBarry Smith     }
206917ab2063SBarry Smith   }
207017ab2063SBarry Smith 
207117ab2063SBarry Smith   /* read in nonzero values */
20720752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr);
207317ab2063SBarry Smith 
207417ab2063SBarry Smith   /* set matrix "i" values */
2075416022c9SBarry Smith   a->i[0] = -shift;
207617ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
2077416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2078416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
207917ab2063SBarry Smith   }
20800452661fSBarry Smith   PetscFree(rowlengths);
208117ab2063SBarry Smith 
20826d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20836d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20843a40ed3dSBarry Smith   PetscFunctionReturn(0);
208517ab2063SBarry Smith }
208617ab2063SBarry Smith 
20875615d1e5SSatish Balay #undef __FUNC__
20885615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
20898f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
20907264ac53SSatish Balay {
20917264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
20927264ac53SSatish Balay 
20933a40ed3dSBarry Smith   PetscFunctionBegin;
2094a8c6a408SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
20957264ac53SSatish Balay 
20967264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
20977264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
2098bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
20993a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2100bcd2baecSBarry Smith   }
21017264ac53SSatish Balay 
21027264ac53SSatish Balay   /* if the a->i are the same */
21038108c231SLois Curfman McInnes   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
21043a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
21057264ac53SSatish Balay   }
21067264ac53SSatish Balay 
21077264ac53SSatish Balay   /* if a->j are the same */
2108bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
21093a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2110bcd2baecSBarry Smith   }
2111bcd2baecSBarry Smith 
2112bcd2baecSBarry Smith   /* if a->a are the same */
211319bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
21143a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
21157264ac53SSatish Balay   }
211677c4ece6SBarry Smith   *flg = PETSC_TRUE;
21173a40ed3dSBarry Smith   PetscFunctionReturn(0);
21187264ac53SSatish Balay 
21197264ac53SSatish Balay }
2120