xref: /petsc/src/mat/impls/aij/seq/aij.c (revision e20fef11c1a1457fc77d0de77d14a2dd6b3370e4)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*e20fef11SSatish Balay static char vcid[] = "$Id: aij.c,v 1.269 1998/05/29 20:37:05 bsmith Exp balay $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5d5d45c9bSBarry Smith /*
63369ce9aSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7d5d45c9bSBarry Smith   matrix storage format.
8d5d45c9bSBarry Smith */
93369ce9aSBarry Smith 
103369ce9aSBarry Smith #include "pinclude/pviewer.h"
113369ce9aSBarry Smith #include "sys.h"
1270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
13f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
14f5eb4b81SSatish Balay #include "src/inline/spops.h"
158d195f9aSBarry Smith #include "src/inline/dot.h"
16f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1717ab2063SBarry Smith 
18a2ce50c7SBarry Smith /*
19a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
20a2ce50c7SBarry Smith */
215615d1e5SSatish Balay #undef __FUNC__
225615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ"
23a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
24a2ce50c7SBarry Smith {
25a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
26a2ce50c7SBarry Smith   int        ierr = 1;
27a2ce50c7SBarry Smith 
283a40ed3dSBarry Smith   PetscFunctionBegin;
29e3372554SBarry Smith   SETERRQ(ierr,0,"Not implemented");
30d64ed03dSBarry Smith #if !defined(USE_PETSC_DEBUG)
31d64ed03dSBarry Smith   PetscFunctionReturn(0);
32d64ed03dSBarry Smith #endif
33a2ce50c7SBarry Smith }
34a2ce50c7SBarry Smith 
35bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3617ab2063SBarry Smith 
375615d1e5SSatish Balay #undef __FUNC__
385615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ"
398f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
406945ee14SBarry Smith                            PetscTruth *done)
4117ab2063SBarry Smith {
42416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
436945ee14SBarry Smith   int        ierr,i,ishift;
4417ab2063SBarry Smith 
453a40ed3dSBarry Smith   PetscFunctionBegin;
4631625ec5SSatish Balay   *m     = A->m;
473a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
486945ee14SBarry Smith   ishift = a->indexshift;
496945ee14SBarry Smith   if (symmetric) {
5031625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
516945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
5231625ec5SSatish Balay     int nz = a->i[a->m];
533b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
5431625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
553b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
563b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
5731625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
586945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
5931625ec5SSatish Balay     int nz = a->i[a->m] + 1;
603b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
6131625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
623b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
633b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
6431625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
656945ee14SBarry Smith   } else {
666945ee14SBarry Smith     *ia = a->i; *ja = a->j;
67a2ce50c7SBarry Smith   }
68a2ce50c7SBarry Smith 
693a40ed3dSBarry Smith   PetscFunctionReturn(0);
70a2744918SBarry Smith }
71a2744918SBarry Smith 
725615d1e5SSatish Balay #undef __FUNC__
735615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ"
748f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
756945ee14SBarry Smith                                PetscTruth *done)
766945ee14SBarry Smith {
776945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
783b2fbd54SBarry Smith   int        ishift = a->indexshift;
796945ee14SBarry Smith 
803a40ed3dSBarry Smith   PetscFunctionBegin;
813a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
823b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
836945ee14SBarry Smith     PetscFree(*ia);
846945ee14SBarry Smith     PetscFree(*ja);
85bcd2baecSBarry Smith   }
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
8717ab2063SBarry Smith }
8817ab2063SBarry Smith 
895615d1e5SSatish Balay #undef __FUNC__
905615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ"
9143a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
923b2fbd54SBarry Smith                            PetscTruth *done)
933b2fbd54SBarry Smith {
943b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
95a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
96a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
973b2fbd54SBarry Smith 
983a40ed3dSBarry Smith   PetscFunctionBegin;
993b2fbd54SBarry Smith   *nn     = A->n;
1003a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1013b2fbd54SBarry Smith   if (symmetric) {
102179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
1033b2fbd54SBarry Smith   } else {
10461d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
1053b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1063b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
107a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
1083b2fbd54SBarry Smith     jj = a->j;
1093b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
1103b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
1113b2fbd54SBarry Smith     }
1123b2fbd54SBarry Smith     cia[0] = oshift;
1133b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
1143b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
1153b2fbd54SBarry Smith     }
1163b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1173b2fbd54SBarry Smith     jj = a->j;
118a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
119a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
120a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1213b2fbd54SBarry Smith         col = *jj++ + ishift;
1223b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1233b2fbd54SBarry Smith       }
1243b2fbd54SBarry Smith     }
1253b2fbd54SBarry Smith     PetscFree(collengths);
1263b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1273b2fbd54SBarry Smith   }
1283b2fbd54SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1303b2fbd54SBarry Smith }
1313b2fbd54SBarry Smith 
1325615d1e5SSatish Balay #undef __FUNC__
1335615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ"
13443a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1353b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1363b2fbd54SBarry Smith {
1373a40ed3dSBarry Smith   PetscFunctionBegin;
1383a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1393b2fbd54SBarry Smith 
1403b2fbd54SBarry Smith   PetscFree(*ia);
1413b2fbd54SBarry Smith   PetscFree(*ja);
1423b2fbd54SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1443b2fbd54SBarry Smith }
1453b2fbd54SBarry Smith 
146227d817aSBarry Smith #define CHUNKSIZE   15
14717ab2063SBarry Smith 
1485615d1e5SSatish Balay #undef __FUNC__
1495615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ"
15061d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
15117ab2063SBarry Smith {
152416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
153416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1544b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
155d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
156416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
15717ab2063SBarry Smith 
1583a40ed3dSBarry Smith   PetscFunctionBegin;
15917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
160416022c9SBarry Smith     row  = im[k];
1613a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
162a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
163a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
1643b2fbd54SBarry Smith #endif
16517ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
16617ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
167416022c9SBarry Smith     low = 0;
16817ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
170a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
171a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
1723b2fbd54SBarry Smith #endif
1734b0e389bSBarry Smith       col = in[l] - shift;
1744b0e389bSBarry Smith       if (roworiented) {
1754b0e389bSBarry Smith         value = *v++;
176bef8e0ddSBarry Smith       } else {
1774b0e389bSBarry Smith         value = v[k + l*m];
1784b0e389bSBarry Smith       }
179416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
180416022c9SBarry Smith       while (high-low > 5) {
181416022c9SBarry Smith         t = (low+high)/2;
182416022c9SBarry Smith         if (rp[t] > col) high = t;
183416022c9SBarry Smith         else             low  = t;
18417ab2063SBarry Smith       }
185416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
18617ab2063SBarry Smith         if (rp[i] > col) break;
18717ab2063SBarry Smith         if (rp[i] == col) {
188416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
18917ab2063SBarry Smith           else                  ap[i] = value;
19017ab2063SBarry Smith           goto noinsert;
19117ab2063SBarry Smith         }
19217ab2063SBarry Smith       }
193c2653b3dSLois Curfman McInnes       if (nonew == 1) goto noinsert;
194a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
19517ab2063SBarry Smith       if (nrow >= rmax) {
19617ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
197416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
19817ab2063SBarry Smith         Scalar *new_a;
19917ab2063SBarry Smith 
200a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
20196854ed6SLois Curfman McInnes 
20217ab2063SBarry Smith         /* malloc new storage space */
203416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
2040452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
20517ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
20617ab2063SBarry Smith         new_i   = new_j + new_nz;
20717ab2063SBarry Smith 
20817ab2063SBarry Smith         /* copy over old data into new slots */
20917ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
210416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
211416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
212416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
213416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
21417ab2063SBarry Smith                                                            len*sizeof(int));
215416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
216416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
21717ab2063SBarry Smith                                                            len*sizeof(Scalar));
21817ab2063SBarry Smith         /* free up old matrix storage */
2190452661fSBarry Smith         PetscFree(a->a);
2200452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
221416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
222416022c9SBarry Smith         a->singlemalloc = 1;
22317ab2063SBarry Smith 
22417ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
225416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
226416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
227416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
228b810aeb4SBarry Smith         a->reallocs++;
22917ab2063SBarry Smith       }
230416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
231416022c9SBarry Smith       /* shift up all the later entries in this row */
232416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
23317ab2063SBarry Smith         rp[ii+1] = rp[ii];
23417ab2063SBarry Smith         ap[ii+1] = ap[ii];
23517ab2063SBarry Smith       }
23617ab2063SBarry Smith       rp[i] = col;
23717ab2063SBarry Smith       ap[i] = value;
23817ab2063SBarry Smith       noinsert:;
239416022c9SBarry Smith       low = i + 1;
24017ab2063SBarry Smith     }
24117ab2063SBarry Smith     ailen[row] = nrow;
24217ab2063SBarry Smith   }
2433a40ed3dSBarry Smith   PetscFunctionReturn(0);
24417ab2063SBarry Smith }
24517ab2063SBarry Smith 
2465615d1e5SSatish Balay #undef __FUNC__
2475615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ"
2488f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2497eb43aa7SLois Curfman McInnes {
2507eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
251b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2527eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2537eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2547eb43aa7SLois Curfman McInnes 
2553a40ed3dSBarry Smith   PetscFunctionBegin;
2567eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2577eb43aa7SLois Curfman McInnes     row  = im[k];
258a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
259a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
2607eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2617eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2627eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
263a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
264a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
2657eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2667eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2677eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2687eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2697eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2707eb43aa7SLois Curfman McInnes         else             low  = t;
2717eb43aa7SLois Curfman McInnes       }
2727eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2737eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2747eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
275b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2767eb43aa7SLois Curfman McInnes           goto finished;
2777eb43aa7SLois Curfman McInnes         }
2787eb43aa7SLois Curfman McInnes       }
279b49de8d1SLois Curfman McInnes       *v++ = zero;
2807eb43aa7SLois Curfman McInnes       finished:;
2817eb43aa7SLois Curfman McInnes     }
2827eb43aa7SLois Curfman McInnes   }
2833a40ed3dSBarry Smith   PetscFunctionReturn(0);
2847eb43aa7SLois Curfman McInnes }
2857eb43aa7SLois Curfman McInnes 
28617ab2063SBarry Smith 
2875615d1e5SSatish Balay #undef __FUNC__
2885615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary"
289480ef9eaSBarry Smith int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
29017ab2063SBarry Smith {
291416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
292416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
29317ab2063SBarry Smith 
2943a40ed3dSBarry Smith   PetscFunctionBegin;
29590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2960452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
297416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
298416022c9SBarry Smith   col_lens[1] = a->m;
299416022c9SBarry Smith   col_lens[2] = a->n;
300416022c9SBarry Smith   col_lens[3] = a->nz;
301416022c9SBarry Smith 
302416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
303416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
304416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
30517ab2063SBarry Smith   }
3060752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3070452661fSBarry Smith   PetscFree(col_lens);
308416022c9SBarry Smith 
309416022c9SBarry Smith   /* store column indices (zero start index) */
310416022c9SBarry Smith   if (a->indexshift) {
311416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
31217ab2063SBarry Smith   }
3130752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0); CHKERRQ(ierr);
314416022c9SBarry Smith   if (a->indexshift) {
315416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
31617ab2063SBarry Smith   }
317416022c9SBarry Smith 
318416022c9SBarry Smith   /* store nonzero values */
3190752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
32117ab2063SBarry Smith }
322416022c9SBarry Smith 
3235615d1e5SSatish Balay #undef __FUNC__
3245615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII"
325480ef9eaSBarry Smith int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
326416022c9SBarry Smith {
327416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
328496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
32917ab2063SBarry Smith   FILE        *fd;
33017ab2063SBarry Smith   char        *outputname;
33117ab2063SBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
33390ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
334416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
33590ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
336a93ec695SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
3373a40ed3dSBarry Smith     PetscFunctionReturn(0);
3383a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
339496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
340496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
341496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
34295e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
34395e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
3443a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
345d00d2cf4SBarry Smith     int nofinalvalue = 0;
346d00d2cf4SBarry Smith     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) {
347d00d2cf4SBarry Smith       nofinalvalue = 1;
348d00d2cf4SBarry Smith     }
349416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3504e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
351d00d2cf4SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue);
35217ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
35317ab2063SBarry Smith 
35417ab2063SBarry Smith     for (i=0; i<m; i++) {
355416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
3563a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
357*e20fef11SSatish 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)
372*e20fef11SSatish Balay         if (PetscImaginary(a->a[j]) > 0.0 && PetscReal(a->a[j]) != 0.0)
373*e20fef11SSatish Balay           fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
374*e20fef11SSatish Balay         else if (PetscImaginary(a->a[j]) < 0.0 && PetscReal(a->a[j]) != 0.0)
375*e20fef11SSatish Balay           fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j]));
376*e20fef11SSatish Balay         else if (PetscReal(a->a[j]) != 0.0)
377*e20fef11SSatish 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)
393*e20fef11SSatish 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)
423*e20fef11SSatish Balay           if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0)
424*e20fef11SSatish 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)
437*e20fef11SSatish Balay         if (PetscImaginary(a->a[j]) > 0.0) {
438*e20fef11SSatish Balay           fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
439*e20fef11SSatish Balay         } else if (PetscImaginary(a->a[j]) < 0.0) {
440*e20fef11SSatish Balay           fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j]));
4413a40ed3dSBarry Smith         } else {
442*e20fef11SSatish 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)
481*e20fef11SSatish 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)
503*e20fef11SSatish 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   }
6703a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
671e1311b90SBarry Smith   PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
67217ab2063SBarry Smith #endif
6730452661fSBarry Smith   PetscFree(a->a);
6740452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
6750452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
6760452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
6770452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
6780452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
67976dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
68082bf6240SBarry Smith   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
6810452661fSBarry Smith   PetscFree(a);
682eed86810SBarry Smith 
683f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
684f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
6853a40ed3dSBarry Smith   PetscFunctionReturn(0);
68617ab2063SBarry Smith }
68717ab2063SBarry Smith 
6885615d1e5SSatish Balay #undef __FUNC__
6895615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
6908f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A)
69117ab2063SBarry Smith {
6923a40ed3dSBarry Smith   PetscFunctionBegin;
6933a40ed3dSBarry Smith   PetscFunctionReturn(0);
69417ab2063SBarry Smith }
69517ab2063SBarry Smith 
6965615d1e5SSatish Balay #undef __FUNC__
6975615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
6988f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op)
69917ab2063SBarry Smith {
700416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7013a40ed3dSBarry Smith 
7023a40ed3dSBarry Smith   PetscFunctionBegin;
7036d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
7046d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
7056d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
706219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
7076d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
708c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
70996854ed6SLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
7106d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
7116d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
712219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7136d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7146d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
71590f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
716b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES||
717b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
718981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n");
7193a40ed3dSBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS) {
7203a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7213a40ed3dSBarry Smith   } else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
7226d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
7236d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
7246d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
7256d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
7263a40ed3dSBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7273a40ed3dSBarry Smith   PetscFunctionReturn(0);
72817ab2063SBarry Smith }
72917ab2063SBarry Smith 
7305615d1e5SSatish Balay #undef __FUNC__
7315615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
7328f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
73317ab2063SBarry Smith {
734416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7353a40ed3dSBarry Smith   int        i,j, n,shift = a->indexshift,ierr;
73617ab2063SBarry Smith   Scalar     *x, zero = 0.0;
73717ab2063SBarry Smith 
7383a40ed3dSBarry Smith   PetscFunctionBegin;
7393a40ed3dSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
740e1311b90SBarry Smith   ierr = VecGetArray(v,&x);;CHKERRQ(ierr);
741e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr);
742a8c6a408SBarry Smith   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
743416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
744416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
745416022c9SBarry Smith       if (a->j[j]+shift == i) {
746416022c9SBarry Smith         x[i] = a->a[j];
74717ab2063SBarry Smith         break;
74817ab2063SBarry Smith       }
74917ab2063SBarry Smith     }
75017ab2063SBarry Smith   }
751e1311b90SBarry Smith   ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr);
7523a40ed3dSBarry Smith   PetscFunctionReturn(0);
75317ab2063SBarry Smith }
75417ab2063SBarry Smith 
75517ab2063SBarry Smith /* -------------------------------------------------------*/
75617ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
75717ab2063SBarry Smith /* -------------------------------------------------------*/
7585615d1e5SSatish Balay #undef __FUNC__
7595615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
76044cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
76117ab2063SBarry Smith {
762416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
76317ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
764e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx, shift = a->indexshift;
76517ab2063SBarry Smith 
7663a40ed3dSBarry Smith   PetscFunctionBegin;
767e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
768e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
769cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
77017ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
77117ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
772416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
773416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
774416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
77517ab2063SBarry Smith     alpha = x[i];
77617ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
77717ab2063SBarry Smith   }
778416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
779e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
780e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7813a40ed3dSBarry Smith   PetscFunctionReturn(0);
78217ab2063SBarry Smith }
783d5d45c9bSBarry Smith 
7845615d1e5SSatish Balay #undef __FUNC__
7855615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
78644cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
78717ab2063SBarry Smith {
788416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
78917ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
790e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx,shift = a->indexshift;
79117ab2063SBarry Smith 
7923a40ed3dSBarry Smith   PetscFunctionBegin;
793e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
794e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
79517ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
79617ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
79717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
798416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
799416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
800416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
80117ab2063SBarry Smith     alpha = x[i];
80217ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
80317ab2063SBarry Smith   }
80490f02eecSBarry Smith   PLogFlops(2*a->nz);
805e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
806e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8073a40ed3dSBarry Smith   PetscFunctionReturn(0);
80817ab2063SBarry Smith }
80917ab2063SBarry Smith 
8105615d1e5SSatish Balay #undef __FUNC__
8115615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
81244cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
81317ab2063SBarry Smith {
814416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
81517ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
816e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
817e36a17ebSSatish Balay #if !defined(USE_FORTRAN_KERNELS)
818e36a17ebSSatish Balay   int        n, i, jrow,j;
819e36a17ebSSatish Balay #endif
82017ab2063SBarry Smith 
821fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT)
822fee21e36SBarry Smith #pragma disjoint(*x,*y,*v)
823fee21e36SBarry Smith #endif
824fee21e36SBarry Smith 
8253a40ed3dSBarry Smith   PetscFunctionBegin;
826e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
827e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
82817ab2063SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
829416022c9SBarry Smith   idx  = a->j;
830d64ed03dSBarry Smith   v    = a->a;
831416022c9SBarry Smith   ii   = a->i;
8328d195f9aSBarry Smith #if defined(USE_FORTRAN_KERNELS)
83329b5ca8aSSatish Balay   fortranmultaij_(&m,x,ii,idx+shift,v+shift,y);
8348d195f9aSBarry Smith #else
835d64ed03dSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
836d64ed03dSBarry Smith   idx  += shift;
83717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8389ea0dfa2SSatish Balay     jrow = ii[i];
8399ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
84017ab2063SBarry Smith     sum  = 0.0;
8419ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8429ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8439ea0dfa2SSatish Balay      }
84417ab2063SBarry Smith     y[i] = sum;
84517ab2063SBarry Smith   }
8468d195f9aSBarry Smith #endif
847416022c9SBarry Smith   PLogFlops(2*a->nz - m);
848e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
849e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
8503a40ed3dSBarry Smith   PetscFunctionReturn(0);
85117ab2063SBarry Smith }
85217ab2063SBarry Smith 
8535615d1e5SSatish Balay #undef __FUNC__
8545615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
85544cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
85617ab2063SBarry Smith {
857416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
85817ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
859e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
860e36a17ebSSatish Balay #if !defined(USE_FORTRAN_KERNELS)
861e36a17ebSSatish Balay   int        n,i,jrow,j;
862e36a17ebSSatish Balay #endif
8639ea0dfa2SSatish Balay 
8643a40ed3dSBarry Smith   PetscFunctionBegin;
865e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
866e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
867e1311b90SBarry Smith   ierr = VecGetArray(zz,&z); CHKERRQ(ierr);
86817ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
869cddf8d76SBarry Smith   idx  = a->j;
870d64ed03dSBarry Smith   v    = a->a;
871cddf8d76SBarry Smith   ii   = a->i;
87202ab625aSSatish Balay #if defined(USE_FORTRAN_KERNELS)
87329b5ca8aSSatish Balay   fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
87402ab625aSSatish Balay #else
875d64ed03dSBarry Smith   v   += shift; /* shift for Fortran start by 1 indexing */
876d64ed03dSBarry Smith   idx += shift;
87717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8789ea0dfa2SSatish Balay     jrow = ii[i];
8799ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
88017ab2063SBarry Smith     sum  = y[i];
8819ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8829ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8839ea0dfa2SSatish Balay      }
88417ab2063SBarry Smith     z[i] = sum;
88517ab2063SBarry Smith   }
88602ab625aSSatish Balay #endif
887416022c9SBarry Smith   PLogFlops(2*a->nz);
888e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
889e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
890e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr);
8913a40ed3dSBarry Smith   PetscFunctionReturn(0);
89217ab2063SBarry Smith }
89317ab2063SBarry Smith 
89417ab2063SBarry Smith /*
89517ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
89617ab2063SBarry Smith */
89717ab2063SBarry Smith 
8985615d1e5SSatish Balay #undef __FUNC__
8995615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
900416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
90117ab2063SBarry Smith {
902416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
903416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
90417ab2063SBarry Smith 
9053a40ed3dSBarry Smith   PetscFunctionBegin;
9060452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
907416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
908416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
909416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
910416022c9SBarry Smith       if (a->j[j]+shift == i) {
91117ab2063SBarry Smith         diag[i] = j - shift;
91217ab2063SBarry Smith         break;
91317ab2063SBarry Smith       }
91417ab2063SBarry Smith     }
91517ab2063SBarry Smith   }
916416022c9SBarry Smith   a->diag = diag;
9173a40ed3dSBarry Smith   PetscFunctionReturn(0);
91817ab2063SBarry Smith }
91917ab2063SBarry Smith 
9205615d1e5SSatish Balay #undef __FUNC__
9215615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
92244cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
92317ab2063SBarry Smith                            double fshift,int its,Vec xx)
92417ab2063SBarry Smith {
925416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
926416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
927d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
92817ab2063SBarry Smith 
9293a40ed3dSBarry Smith   PetscFunctionBegin;
930e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
931e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
932d00d2cf4SBarry Smith   if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);}
933416022c9SBarry Smith   diag = a->diag;
934416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
93517ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
93617ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
93717ab2063SBarry Smith     bs = b + shift;
93817ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
939416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
940416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
941416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
942416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
94317ab2063SBarry Smith         sum  = b[i]*d/omega;
94417ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
94517ab2063SBarry Smith         x[i] = sum;
94617ab2063SBarry Smith     }
9473a40ed3dSBarry Smith     PetscFunctionReturn(0);
94817ab2063SBarry Smith   }
94917ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
950a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done");
9513a40ed3dSBarry Smith   } else if (flag & SOR_EISENSTAT) {
95217ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
95317ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
95417ab2063SBarry Smith 
95517ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
95617ab2063SBarry Smith 
95717ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
95817ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
95917ab2063SBarry Smith     is the relaxation factor.
96017ab2063SBarry Smith     */
9610452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
96217ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
96317ab2063SBarry Smith 
96417ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
96517ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
966416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
967416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
968416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
969416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
97017ab2063SBarry Smith       sum  = b[i];
97117ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
97217ab2063SBarry Smith       x[i] = omega*(sum/d);
97317ab2063SBarry Smith     }
97417ab2063SBarry Smith 
97517ab2063SBarry Smith     /*  t = b - (2*E - D)x */
976416022c9SBarry Smith     v = a->a;
97717ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
97817ab2063SBarry Smith 
97917ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
980416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
981416022c9SBarry Smith     diag = a->diag;
98217ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
983416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
984416022c9SBarry Smith       n    = diag[i] - a->i[i];
985416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
986416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
98717ab2063SBarry Smith       sum  = t[i];
98817ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
98917ab2063SBarry Smith       t[i] = omega*(sum/d);
99017ab2063SBarry Smith     }
99117ab2063SBarry Smith 
99217ab2063SBarry Smith     /*  x = x + t */
99317ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
9940452661fSBarry Smith     PetscFree(t);
9953a40ed3dSBarry Smith     PetscFunctionReturn(0);
99617ab2063SBarry Smith   }
99717ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
99817ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
99917ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1000416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1001416022c9SBarry Smith         n    = diag[i] - a->i[i];
1002416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1003416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
100417ab2063SBarry Smith         sum  = b[i];
100517ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
100617ab2063SBarry Smith         x[i] = omega*(sum/d);
100717ab2063SBarry Smith       }
100817ab2063SBarry Smith       xb = x;
10093a40ed3dSBarry Smith     } else xb = b;
101017ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
101117ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
101217ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1013416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
101417ab2063SBarry Smith       }
101517ab2063SBarry Smith     }
101617ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
101717ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1018416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1019416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
1020416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1021416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
102217ab2063SBarry Smith         sum  = xb[i];
102317ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
102417ab2063SBarry Smith         x[i] = omega*(sum/d);
102517ab2063SBarry Smith       }
102617ab2063SBarry Smith     }
102717ab2063SBarry Smith     its--;
102817ab2063SBarry Smith   }
102917ab2063SBarry Smith   while (its--) {
103017ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
103117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1032416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1033416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1034416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1035416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
103617ab2063SBarry Smith         sum  = b[i];
103717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10387e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
103917ab2063SBarry Smith       }
104017ab2063SBarry Smith     }
104117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
104217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1043416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1044416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1045416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1046416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
104717ab2063SBarry Smith         sum  = b[i];
104817ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10497e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
105017ab2063SBarry Smith       }
105117ab2063SBarry Smith     }
105217ab2063SBarry Smith   }
10533a40ed3dSBarry Smith   PetscFunctionReturn(0);
105417ab2063SBarry Smith }
105517ab2063SBarry Smith 
10565615d1e5SSatish Balay #undef __FUNC__
10575615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10588f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
105917ab2063SBarry Smith {
1060416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
10614e220ebcSLois Curfman McInnes 
10623a40ed3dSBarry Smith   PetscFunctionBegin;
10634e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
10644e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
10654e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
10664e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
10674e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10684e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
10694e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
10704e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
10714e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
10724e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
10734e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
10744e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
10754e220ebcSLois Curfman McInnes   info->memory         = A->mem;
10764e220ebcSLois Curfman McInnes   if (A->factor) {
10774e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10784e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
10794e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
10804e220ebcSLois Curfman McInnes   } else {
10814e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
10824e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
10834e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
10844e220ebcSLois Curfman McInnes   }
10853a40ed3dSBarry Smith   PetscFunctionReturn(0);
108617ab2063SBarry Smith }
108717ab2063SBarry Smith 
108817ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
108917ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
109017ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
109117ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
109217ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
109317ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
109417ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
109517ab2063SBarry Smith 
10965615d1e5SSatish Balay #undef __FUNC__
10975615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
10988f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
109917ab2063SBarry Smith {
1100416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1101416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
110217ab2063SBarry Smith 
11033a40ed3dSBarry Smith   PetscFunctionBegin;
110477c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
110517ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
110617ab2063SBarry Smith   if (diag) {
110717ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1108a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1109416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1110416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1111416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1112416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
11133a40ed3dSBarry Smith       } else {
1114d64ed03dSBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr);
111517ab2063SBarry Smith       }
111617ab2063SBarry Smith     }
11173a40ed3dSBarry Smith   } else {
111817ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1119a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1120416022c9SBarry Smith       a->ilen[rows[i]] = 0;
112117ab2063SBarry Smith     }
112217ab2063SBarry Smith   }
112317ab2063SBarry Smith   ISRestoreIndices(is,&rows);
112443a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11253a40ed3dSBarry Smith   PetscFunctionReturn(0);
112617ab2063SBarry Smith }
112717ab2063SBarry Smith 
11285615d1e5SSatish Balay #undef __FUNC__
11295615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11308f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
113117ab2063SBarry Smith {
1132416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11333a40ed3dSBarry Smith 
11343a40ed3dSBarry Smith   PetscFunctionBegin;
11350752156aSBarry Smith   if (m) *m = a->m;
11360752156aSBarry Smith   if (n) *n = a->n;
11373a40ed3dSBarry Smith   PetscFunctionReturn(0);
113817ab2063SBarry Smith }
113917ab2063SBarry Smith 
11405615d1e5SSatish Balay #undef __FUNC__
11415615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11428f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
114317ab2063SBarry Smith {
1144416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11453a40ed3dSBarry Smith 
11463a40ed3dSBarry Smith   PetscFunctionBegin;
1147416022c9SBarry Smith   *m = 0; *n = a->m;
11483a40ed3dSBarry Smith   PetscFunctionReturn(0);
114917ab2063SBarry Smith }
1150026e39d0SSatish Balay 
11515615d1e5SSatish Balay #undef __FUNC__
11525615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11534e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
115417ab2063SBarry Smith {
1155416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1156c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
115717ab2063SBarry Smith 
11583a40ed3dSBarry Smith   PetscFunctionBegin;
1159a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
116017ab2063SBarry Smith 
1161416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1162416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
116317ab2063SBarry Smith   if (idx) {
1164416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
11654e093b46SBarry Smith     if (*nz && shift) {
11660452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
116717ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
11684e093b46SBarry Smith     } else if (*nz) {
11694e093b46SBarry Smith       *idx = itmp;
117017ab2063SBarry Smith     }
117117ab2063SBarry Smith     else *idx = 0;
117217ab2063SBarry Smith   }
11733a40ed3dSBarry Smith   PetscFunctionReturn(0);
117417ab2063SBarry Smith }
117517ab2063SBarry Smith 
11765615d1e5SSatish Balay #undef __FUNC__
11775615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
11784e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
117917ab2063SBarry Smith {
11804e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11813a40ed3dSBarry Smith 
11823a40ed3dSBarry Smith   PetscFunctionBegin;
11834e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
11843a40ed3dSBarry Smith   PetscFunctionReturn(0);
118517ab2063SBarry Smith }
118617ab2063SBarry Smith 
11875615d1e5SSatish Balay #undef __FUNC__
11885615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
11898f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
119017ab2063SBarry Smith {
1191416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1192416022c9SBarry Smith   Scalar     *v = a->a;
119317ab2063SBarry Smith   double     sum = 0.0;
1194416022c9SBarry Smith   int        i, j,shift = a->indexshift;
119517ab2063SBarry Smith 
11963a40ed3dSBarry Smith   PetscFunctionBegin;
119717ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1198416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
11993a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
1200*e20fef11SSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
120117ab2063SBarry Smith #else
120217ab2063SBarry Smith       sum += (*v)*(*v); v++;
120317ab2063SBarry Smith #endif
120417ab2063SBarry Smith     }
120517ab2063SBarry Smith     *norm = sqrt(sum);
12063a40ed3dSBarry Smith   } else if (type == NORM_1) {
120717ab2063SBarry Smith     double *tmp;
1208416022c9SBarry Smith     int    *jj = a->j;
120966963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1210cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
121117ab2063SBarry Smith     *norm = 0.0;
1212416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1213a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
121417ab2063SBarry Smith     }
1215416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
121617ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
121717ab2063SBarry Smith     }
12180452661fSBarry Smith     PetscFree(tmp);
12193a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
122017ab2063SBarry Smith     *norm = 0.0;
1221416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1222416022c9SBarry Smith       v = a->a + a->i[j] + shift;
122317ab2063SBarry Smith       sum = 0.0;
1224416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1225cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
122617ab2063SBarry Smith       }
122717ab2063SBarry Smith       if (sum > *norm) *norm = sum;
122817ab2063SBarry Smith     }
12293a40ed3dSBarry Smith   } else {
1230a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
123117ab2063SBarry Smith   }
12323a40ed3dSBarry Smith   PetscFunctionReturn(0);
123317ab2063SBarry Smith }
123417ab2063SBarry Smith 
12355615d1e5SSatish Balay #undef __FUNC__
12365615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12378f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
123817ab2063SBarry Smith {
1239416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1240416022c9SBarry Smith   Mat        C;
1241416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1242416022c9SBarry Smith   int        shift = a->indexshift;
1243d5d45c9bSBarry Smith   Scalar     *array = a->a;
124417ab2063SBarry Smith 
12453a40ed3dSBarry Smith   PetscFunctionBegin;
1246a8c6a408SBarry Smith   if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place");
12470452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1248cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
124917ab2063SBarry Smith   if (shift) {
125017ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
125117ab2063SBarry Smith   }
125217ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1253416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12540452661fSBarry Smith   PetscFree(col);
125517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
125617ab2063SBarry Smith     len = ai[i+1]-ai[i];
1257416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
125817ab2063SBarry Smith     array += len; aj += len;
125917ab2063SBarry Smith   }
126017ab2063SBarry Smith   if (shift) {
126117ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
126217ab2063SBarry Smith   }
126317ab2063SBarry Smith 
12646d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12656d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
126617ab2063SBarry Smith 
12673638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1268416022c9SBarry Smith     *B = C;
126917ab2063SBarry Smith   } else {
1270f830108cSBarry Smith     PetscOps       *Abops;
1271f830108cSBarry Smith     struct _MatOps *Aops;
1272f830108cSBarry Smith 
1273416022c9SBarry Smith     /* This isn't really an in-place transpose */
12740452661fSBarry Smith     PetscFree(a->a);
12750452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
12760452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
12770452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
12780452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
12790452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
12801073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
12810452661fSBarry Smith     PetscFree(a);
1282f830108cSBarry Smith 
1283f830108cSBarry Smith     /*
1284f830108cSBarry Smith       This is horrible, horrible code. We need to keep the
12858f0f457cSSatish Balay       the bops and ops Structures, copy everything from C
12868f0f457cSSatish Balay       including the function pointers..
1287f830108cSBarry Smith     */
12888f0f457cSSatish Balay     PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps));
12898f0f457cSSatish Balay     PetscMemcpy(A->bops,C->bops,sizeof(PetscOps));
1290f830108cSBarry Smith     Abops    = A->bops;
1291f830108cSBarry Smith     Aops     = A->ops;
1292f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
1293f830108cSBarry Smith     A->bops  = Abops;
1294f830108cSBarry Smith     A->ops   = Aops;
129527a8da17SBarry Smith     A->qlist = 0;
1296f830108cSBarry Smith 
12970452661fSBarry Smith     PetscHeaderDestroy(C);
129817ab2063SBarry Smith   }
12993a40ed3dSBarry Smith   PetscFunctionReturn(0);
130017ab2063SBarry Smith }
130117ab2063SBarry Smith 
13025615d1e5SSatish Balay #undef __FUNC__
13035615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
13048f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
130517ab2063SBarry Smith {
1306416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
130717ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1308e1311b90SBarry Smith   int        ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
130917ab2063SBarry Smith 
13103a40ed3dSBarry Smith   PetscFunctionBegin;
131117ab2063SBarry Smith   if (ll) {
13123ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
13133ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
1314e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr);
1315a8c6a408SBarry Smith     if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
1316e1311b90SBarry Smith     ierr = VecGetArray(ll,&l); CHKERRQ(ierr);
1317416022c9SBarry Smith     v = a->a;
131817ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
131917ab2063SBarry Smith       x = l[i];
1320416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
132117ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
132217ab2063SBarry Smith     }
1323e1311b90SBarry Smith     ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr);
132444cd7ae7SLois Curfman McInnes     PLogFlops(nz);
132517ab2063SBarry Smith   }
132617ab2063SBarry Smith   if (rr) {
1327e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr);
1328a8c6a408SBarry Smith     if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
1329e1311b90SBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1330416022c9SBarry Smith     v = a->a; jj = a->j;
133117ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
133217ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
133317ab2063SBarry Smith     }
1334e1311b90SBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
133544cd7ae7SLois Curfman McInnes     PLogFlops(nz);
133617ab2063SBarry Smith   }
13373a40ed3dSBarry Smith   PetscFunctionReturn(0);
133817ab2063SBarry Smith }
133917ab2063SBarry Smith 
13405615d1e5SSatish Balay #undef __FUNC__
13415615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
13426a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B)
134317ab2063SBarry Smith {
1344db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1345d5db1b72SBarry Smith   int          *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
134699141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1347a2744918SBarry Smith   register int sum,lensi;
134899141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
134999141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
135099141d43SSatish Balay   Scalar       *a_new,*mat_a;
1351416022c9SBarry Smith   Mat          C;
1352fee21e36SBarry Smith   PetscTruth   stride;
135317ab2063SBarry Smith 
13543a40ed3dSBarry Smith   PetscFunctionBegin;
1355d64ed03dSBarry Smith   ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr);
1356a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
1357d64ed03dSBarry Smith   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
1358a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
135999141d43SSatish Balay 
136017ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
136117ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
136217ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
136317ab2063SBarry Smith 
1364fee21e36SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr);
1365fee21e36SBarry Smith   ierr = ISStride(iscol,&stride); CHKERRQ(ierr);
1366fee21e36SBarry Smith   if (stride && step == 1) {
136702834360SBarry Smith     /* special case of contiguous rows */
136857faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
136902834360SBarry Smith     starts = lens + ncols;
137002834360SBarry Smith     /* loop over new rows determining lens and starting points */
137102834360SBarry Smith     for (i=0; i<nrows; i++) {
1372a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1373a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
137402834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1375d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
137602834360SBarry Smith           starts[i] = k;
137702834360SBarry Smith           break;
137802834360SBarry Smith 	}
137902834360SBarry Smith       }
1380a2744918SBarry Smith       sum = 0;
138102834360SBarry Smith       while (k < kend) {
1382d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1383a2744918SBarry Smith         sum++;
138402834360SBarry Smith       }
1385a2744918SBarry Smith       lens[i] = sum;
138602834360SBarry Smith     }
138702834360SBarry Smith     /* create submatrix */
1388cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
138908480c60SBarry Smith       int n_cols,n_rows;
139008480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1391a8c6a408SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1392d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
139308480c60SBarry Smith       C = *B;
13943a40ed3dSBarry Smith     } else {
139502834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
139608480c60SBarry Smith     }
1397db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1398db02288aSLois Curfman McInnes 
139902834360SBarry Smith     /* loop over rows inserting into submatrix */
1400db02288aSLois Curfman McInnes     a_new    = c->a;
1401db02288aSLois Curfman McInnes     j_new    = c->j;
1402db02288aSLois Curfman McInnes     i_new    = c->i;
1403db02288aSLois Curfman McInnes     i_new[0] = -shift;
140402834360SBarry Smith     for (i=0; i<nrows; i++) {
1405a2744918SBarry Smith       ii    = starts[i];
1406a2744918SBarry Smith       lensi = lens[i];
1407a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1408a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
140902834360SBarry Smith       }
1410a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1411a2744918SBarry Smith       a_new      += lensi;
1412a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1413a2744918SBarry Smith       c->ilen[i]  = lensi;
141402834360SBarry Smith     }
14150452661fSBarry Smith     PetscFree(lens);
14163a40ed3dSBarry Smith   } else {
141702834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
14180452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
141902834360SBarry Smith     ssmap = smap + shift;
142099141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1421cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
142217ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
142302834360SBarry Smith     /* determine lens of each row */
142402834360SBarry Smith     for (i=0; i<nrows; i++) {
1425d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
142602834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
142702834360SBarry Smith       lens[i] = 0;
142802834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1429d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
143002834360SBarry Smith           lens[i]++;
143102834360SBarry Smith         }
143202834360SBarry Smith       }
143302834360SBarry Smith     }
143417ab2063SBarry Smith     /* Create and fill new matrix */
1435a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
143699141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
1437a8c6a408SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
143899141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1439a8c6a408SBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
144099141d43SSatish Balay       }
144199141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
144208480c60SBarry Smith       C = *B;
14433a40ed3dSBarry Smith     } else {
144402834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
144508480c60SBarry Smith     }
144699141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
144717ab2063SBarry Smith     for (i=0; i<nrows; i++) {
144899141d43SSatish Balay       row    = irow[i];
144999141d43SSatish Balay       kstart = ai[row]+shift;
145099141d43SSatish Balay       kend   = kstart + a->ilen[row];
145199141d43SSatish Balay       mat_i  = c->i[i]+shift;
145299141d43SSatish Balay       mat_j  = c->j + mat_i;
145399141d43SSatish Balay       mat_a  = c->a + mat_i;
145499141d43SSatish Balay       mat_ilen = c->ilen + i;
145517ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
145699141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
145799141d43SSatish Balay           *mat_j++ = tcol - (!shift);
145899141d43SSatish Balay           *mat_a++ = a->a[k];
145999141d43SSatish Balay           (*mat_ilen)++;
146099141d43SSatish Balay 
146117ab2063SBarry Smith         }
146217ab2063SBarry Smith       }
146317ab2063SBarry Smith     }
146402834360SBarry Smith     /* Free work space */
146502834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
146699141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
146702834360SBarry Smith   }
14686d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
14696d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
147017ab2063SBarry Smith 
147117ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1472416022c9SBarry Smith   *B = C;
14733a40ed3dSBarry Smith   PetscFunctionReturn(0);
147417ab2063SBarry Smith }
147517ab2063SBarry Smith 
1476a871dcd8SBarry Smith /*
147763b91edcSBarry Smith      note: This can only work for identity for row and col. It would
147863b91edcSBarry Smith    be good to check this and otherwise generate an error.
1479a871dcd8SBarry Smith */
14805615d1e5SSatish Balay #undef __FUNC__
14815615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
14828f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1483a871dcd8SBarry Smith {
148463b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
148508480c60SBarry Smith   int        ierr;
148663b91edcSBarry Smith   Mat        outA;
148763b91edcSBarry Smith 
14883a40ed3dSBarry Smith   PetscFunctionBegin;
1489a8c6a408SBarry Smith   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
1490a871dcd8SBarry Smith 
149163b91edcSBarry Smith   outA          = inA;
149263b91edcSBarry Smith   inA->factor   = FACTOR_LU;
149363b91edcSBarry Smith   a->row        = row;
149463b91edcSBarry Smith   a->col        = col;
149563b91edcSBarry Smith 
1496f0ec6fceSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1497f0ec6fceSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
14981e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1499f0ec6fceSSatish Balay 
150094a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
15010452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work);
150294a9d846SBarry Smith   }
150363b91edcSBarry Smith 
150408480c60SBarry Smith   if (!a->diag) {
150508480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
150663b91edcSBarry Smith   }
150763b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
15083a40ed3dSBarry Smith   PetscFunctionReturn(0);
1509a871dcd8SBarry Smith }
1510a871dcd8SBarry Smith 
151174cdf0dfSBarry Smith #include "pinclude/blaslapack.h"
15125615d1e5SSatish Balay #undef __FUNC__
15135615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
15148f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1515f0b747eeSBarry Smith {
1516f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1517f0b747eeSBarry Smith   int        one = 1;
15183a40ed3dSBarry Smith 
15193a40ed3dSBarry Smith   PetscFunctionBegin;
1520f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1521f0b747eeSBarry Smith   PLogFlops(a->nz);
15223a40ed3dSBarry Smith   PetscFunctionReturn(0);
1523f0b747eeSBarry Smith }
1524f0b747eeSBarry Smith 
15255615d1e5SSatish Balay #undef __FUNC__
15265615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
15276a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1528cddf8d76SBarry Smith {
1529cddf8d76SBarry Smith   int ierr,i;
1530cddf8d76SBarry Smith 
15313a40ed3dSBarry Smith   PetscFunctionBegin;
1532cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
15330452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1534cddf8d76SBarry Smith   }
1535cddf8d76SBarry Smith 
1536cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
15376a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1538cddf8d76SBarry Smith   }
15393a40ed3dSBarry Smith   PetscFunctionReturn(0);
1540cddf8d76SBarry Smith }
1541cddf8d76SBarry Smith 
15425615d1e5SSatish Balay #undef __FUNC__
15435615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
15448f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
15455a838052SSatish Balay {
1546f830108cSBarry Smith   PetscFunctionBegin;
15475a838052SSatish Balay   *bs = 1;
15483a40ed3dSBarry Smith   PetscFunctionReturn(0);
15495a838052SSatish Balay }
15505a838052SSatish Balay 
15515615d1e5SSatish Balay #undef __FUNC__
15525615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15538f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15544dcbc457SBarry Smith {
1555e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
155606763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15578a047759SSatish Balay   int        start, end, *ai, *aj;
1558bbd702dbSSatish Balay   BT         table;
1559bbd702dbSSatish Balay 
15603a40ed3dSBarry Smith   PetscFunctionBegin;
15618a047759SSatish Balay   shift = a->indexshift;
1562e4d965acSSatish Balay   m     = a->m;
1563e4d965acSSatish Balay   ai    = a->i;
15648a047759SSatish Balay   aj    = a->j+shift;
15658a047759SSatish Balay 
1566a8c6a408SBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used");
156706763907SSatish Balay 
156806763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1569bbd702dbSSatish Balay   ierr  = BTCreate(m,table); CHKERRQ(ierr);
157006763907SSatish Balay 
1571e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1572b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1573e4d965acSSatish Balay     isz  = 0;
1574bbd702dbSSatish Balay     BTMemzero(m,table);
1575e4d965acSSatish Balay 
1576e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
15774dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
157877c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1579e4d965acSSatish Balay 
1580dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1581e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
1582bbd702dbSSatish Balay       if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];}
15834dcbc457SBarry Smith     }
158406763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
158506763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1586e4d965acSSatish Balay 
158704a348a9SBarry Smith     k = 0;
158804a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
158904a348a9SBarry Smith       n = isz;
159006763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1591e4d965acSSatish Balay         row   = nidx[k];
1592e4d965acSSatish Balay         start = ai[row];
1593e4d965acSSatish Balay         end   = ai[row+1];
159404a348a9SBarry Smith         for ( l = start; l<end ; l++){
15958a047759SSatish Balay           val = aj[l] + shift;
15962a8f2162SSatish Balay           if (!BTLookupSet(table,val)) {nidx[isz++] = val;}
1597e4d965acSSatish Balay         }
1598e4d965acSSatish Balay       }
1599e4d965acSSatish Balay     }
1600029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1601e4d965acSSatish Balay   }
1602bbd702dbSSatish Balay   BTDestroy(table);
160306763907SSatish Balay   PetscFree(nidx);
16043a40ed3dSBarry Smith   PetscFunctionReturn(0);
16054dcbc457SBarry Smith }
160617ab2063SBarry Smith 
16070513a670SBarry Smith /* -------------------------------------------------------------- */
16080513a670SBarry Smith #undef __FUNC__
16090513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
16100513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
16110513a670SBarry Smith {
16120513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
16130513a670SBarry Smith   Scalar     *vwork;
16140513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
16150513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
161656cd22aeSBarry Smith   IS         icolp,irowp;
16170513a670SBarry Smith 
16183a40ed3dSBarry Smith   PetscFunctionBegin;
161956cd22aeSBarry Smith   ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr);
162056cd22aeSBarry Smith   ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr);
162156cd22aeSBarry Smith   ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr);
162256cd22aeSBarry Smith   ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr);
16230513a670SBarry Smith 
16240513a670SBarry Smith   /* determine lengths of permuted rows */
16250513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
16260513a670SBarry Smith   for (i=0; i<m; i++ ) {
16270513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
16280513a670SBarry Smith   }
16290513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
16300513a670SBarry Smith   PetscFree(lens);
16310513a670SBarry Smith 
16320513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
16330513a670SBarry Smith   for (i=0; i<m; i++) {
16340513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16350513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
16360513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
16370513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16380513a670SBarry Smith   }
16390513a670SBarry Smith   PetscFree(cnew);
16400513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
16410513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
164256cd22aeSBarry Smith   ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr);
164356cd22aeSBarry Smith   ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr);
164456cd22aeSBarry Smith   ierr = ISDestroy(irowp); CHKERRQ(ierr);
164556cd22aeSBarry Smith   ierr = ISDestroy(icolp); CHKERRQ(ierr);
16463a40ed3dSBarry Smith   PetscFunctionReturn(0);
16470513a670SBarry Smith }
16480513a670SBarry Smith 
16495615d1e5SSatish Balay #undef __FUNC__
16505615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1651682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1652682d7d0cSBarry Smith {
1653682d7d0cSBarry Smith   static int called = 0;
1654682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1655682d7d0cSBarry Smith 
16563a40ed3dSBarry Smith   PetscFunctionBegin;
16573a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
165876be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
165976be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
166076be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
166176be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodes\n");
166276be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1663682d7d0cSBarry Smith #if defined(HAVE_ESSL)
166476be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1665682d7d0cSBarry Smith #endif
16663a40ed3dSBarry Smith   PetscFunctionReturn(0);
1667682d7d0cSBarry Smith }
16688f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1669a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1670a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1671a93ec695SBarry Smith 
1672682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
167317ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
167417ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1675416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1676416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
167717ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
167817ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
167917ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
168017ab2063SBarry Smith        MatRelax_SeqAIJ,
168117ab2063SBarry Smith        MatTranspose_SeqAIJ,
16827264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1683f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
168417ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
168517ab2063SBarry Smith        MatCompress_SeqAIJ,
168617ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
168717ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
168817ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
168917ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
169094a9d846SBarry Smith        0,0,
16913d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1692cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
16937eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1694682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1695f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
16965a838052SSatish Balay        MatScale_SeqAIJ,0,0,
16976945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
16986945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
16993b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
17003b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
17013b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1702a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1703a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
17040513a670SBarry Smith        MatColoringPatch_SeqAIJ,
17050513a670SBarry Smith        0,
17060513a670SBarry Smith        MatPermute_SeqAIJ};
170717ab2063SBarry Smith 
170817ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
170917ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
171017ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
171117ab2063SBarry Smith 
17125615d1e5SSatish Balay #undef __FUNC__
1713bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ"
1714bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
1715bef8e0ddSBarry Smith {
1716bef8e0ddSBarry Smith   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
1717bef8e0ddSBarry Smith   int        i,nz,n;
1718bef8e0ddSBarry Smith 
1719bef8e0ddSBarry Smith   PetscFunctionBegin;
1720bef8e0ddSBarry Smith   if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing");
1721bef8e0ddSBarry Smith 
1722bef8e0ddSBarry Smith   nz = aij->maxnz;
1723bef8e0ddSBarry Smith   n  = aij->n;
1724bef8e0ddSBarry Smith   for (i=0; i<nz; i++) {
1725bef8e0ddSBarry Smith     aij->j[i] = indices[i];
1726bef8e0ddSBarry Smith   }
1727bef8e0ddSBarry Smith   aij->nz = nz;
1728bef8e0ddSBarry Smith   for ( i=0; i<n; i++ ) {
1729bef8e0ddSBarry Smith     aij->ilen[i] = aij->imax[i];
1730bef8e0ddSBarry Smith   }
1731bef8e0ddSBarry Smith 
1732bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1733bef8e0ddSBarry Smith }
1734bef8e0ddSBarry Smith 
1735bef8e0ddSBarry Smith #undef __FUNC__
1736bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices"
1737bef8e0ddSBarry Smith /*@
1738bef8e0ddSBarry Smith     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
1739bef8e0ddSBarry Smith        in the matrix.
1740bef8e0ddSBarry Smith 
1741bef8e0ddSBarry Smith   Input Parameters:
1742bef8e0ddSBarry Smith +  mat - the SeqAIJ matrix
1743bef8e0ddSBarry Smith -  indices - the column indices
1744bef8e0ddSBarry Smith 
1745bef8e0ddSBarry Smith   Notes:
1746bef8e0ddSBarry Smith     This can be called if you have precomputed the nonzero structure of the
1747bef8e0ddSBarry Smith   matrix and want to provide it to the matrix object to improve the performance
1748bef8e0ddSBarry Smith   of the MatSetValues() operation.
1749bef8e0ddSBarry Smith 
1750bef8e0ddSBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
1751bef8e0ddSBarry Smith   MatCreateSeqAIJ().
1752bef8e0ddSBarry Smith 
1753bef8e0ddSBarry Smith     MUST be called before any calls to MatSetValues();
1754bef8e0ddSBarry Smith 
1755bef8e0ddSBarry Smith @*/
1756bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
1757bef8e0ddSBarry Smith {
1758bef8e0ddSBarry Smith   int ierr,(*f)(Mat,int *);
1759bef8e0ddSBarry Smith 
1760bef8e0ddSBarry Smith   PetscFunctionBegin;
1761bef8e0ddSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1762bef8e0ddSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f);
1763bef8e0ddSBarry Smith          CHKERRQ(ierr);
1764bef8e0ddSBarry Smith   if (f) {
1765bef8e0ddSBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
1766bef8e0ddSBarry Smith   } else {
1767bef8e0ddSBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
1768bef8e0ddSBarry Smith   }
1769bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1770bef8e0ddSBarry Smith }
1771bef8e0ddSBarry Smith 
1772bef8e0ddSBarry Smith #undef __FUNC__
17735615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
177417ab2063SBarry Smith /*@C
1775682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
17760d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
17776e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
17782bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
17792bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
178017ab2063SBarry Smith 
1781db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1782db81eaa0SLois Curfman McInnes 
178317ab2063SBarry Smith    Input Parameters:
1784db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
178517ab2063SBarry Smith .  m - number of rows
178617ab2063SBarry Smith .  n - number of columns
178717ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
1788db81eaa0SLois Curfman McInnes -  nzz - array containing the number of nonzeros in the various rows
17892bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
179017ab2063SBarry Smith 
179117ab2063SBarry Smith    Output Parameter:
1792416022c9SBarry Smith .  A - the matrix
179317ab2063SBarry Smith 
1794b259b22eSLois Curfman McInnes    Notes:
179517ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
179617ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
17970002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
179844cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
179917ab2063SBarry Smith 
180017ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1801a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
18023d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
18036da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
180417ab2063SBarry Smith 
1805682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
18064fca80b9SLois Curfman McInnes    improve numerical efficiency of matrix-vector products and solves. We
1807682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
18086c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
18096c7ebb05SLois Curfman McInnes 
18106c7ebb05SLois Curfman McInnes    Options Database Keys:
1811db81eaa0SLois Curfman McInnes +  -mat_aij_no_inode  - Do not use inodes
1812db81eaa0SLois Curfman McInnes .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
1813db81eaa0SLois Curfman McInnes -  -mat_aij_oneindex - Internally use indexing starting at 1
1814db81eaa0SLois Curfman McInnes         rather than 0.  Note that when calling MatSetValues(),
1815db81eaa0SLois Curfman McInnes         the user still MUST index entries starting at 0!
181617ab2063SBarry Smith 
1817bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices()
181817ab2063SBarry Smith @*/
1819416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
182017ab2063SBarry Smith {
1821416022c9SBarry Smith   Mat        B;
1822416022c9SBarry Smith   Mat_SeqAIJ *b;
18236945ee14SBarry Smith   int        i, len, ierr, flg,size;
18246945ee14SBarry Smith 
18253a40ed3dSBarry Smith   PetscFunctionBegin;
18266945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1827a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1");
1828d5d45c9bSBarry Smith 
1829416022c9SBarry Smith   *A                  = 0;
1830f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView);
1831416022c9SBarry Smith   PLogObjectCreate(B);
18320452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
183344cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1834f830108cSBarry Smith   PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps));
1835e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_SeqAIJ;
1836e1311b90SBarry Smith   B->ops->view             = MatView_SeqAIJ;
1837416022c9SBarry Smith   B->factor           = 0;
1838416022c9SBarry Smith   B->lupivotthreshold = 1.0;
183990f02eecSBarry Smith   B->mapping          = 0;
1840e1311b90SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg); CHKERRQ(ierr);
18417a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
1842e1311b90SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*) &b->ilu_preserve_row_sums);CHKERRQ(ierr);
1843416022c9SBarry Smith   b->row              = 0;
1844416022c9SBarry Smith   b->col              = 0;
184582bf6240SBarry Smith   b->icol             = 0;
1846416022c9SBarry Smith   b->indexshift       = 0;
1847b810aeb4SBarry Smith   b->reallocs         = 0;
184869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
184969957df2SSatish Balay   if (flg) b->indexshift = -1;
185017ab2063SBarry Smith 
185144cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
185244cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
18530452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1854b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
18557b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
18567b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1857416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
185817ab2063SBarry Smith     nz = nz*m;
18593a40ed3dSBarry Smith   } else {
186017ab2063SBarry Smith     nz = 0;
1861416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
186217ab2063SBarry Smith   }
186317ab2063SBarry Smith 
186417ab2063SBarry Smith   /* allocate the matrix space */
1865416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
18660452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1867416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1868cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1869416022c9SBarry Smith   b->i  = b->j + nz;
1870416022c9SBarry Smith   b->singlemalloc = 1;
187117ab2063SBarry Smith 
1872416022c9SBarry Smith   b->i[0] = -b->indexshift;
187317ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1874416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
187517ab2063SBarry Smith   }
187617ab2063SBarry Smith 
1877416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
18780452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1879f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1880416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
188117ab2063SBarry Smith 
1882416022c9SBarry Smith   b->nz               = 0;
1883416022c9SBarry Smith   b->maxnz            = nz;
1884416022c9SBarry Smith   b->sorted           = 0;
1885416022c9SBarry Smith   b->roworiented      = 1;
1886416022c9SBarry Smith   b->nonew            = 0;
1887416022c9SBarry Smith   b->diag             = 0;
1888416022c9SBarry Smith   b->solve_work       = 0;
188971bd300dSLois Curfman McInnes   b->spptr            = 0;
1890754ec7b1SSatish Balay   b->inode.node_count = 0;
1891754ec7b1SSatish Balay   b->inode.size       = 0;
18926c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
18936c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
18944e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
189517ab2063SBarry Smith 
1896416022c9SBarry Smith   *A = B;
18974e220ebcSLois Curfman McInnes 
18984b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
18994b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
190069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
190169957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
19024b14c69eSBarry Smith #endif
190369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
190469957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
190569957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
190669957df2SSatish Balay   if (flg) {
1907a8c6a408SBarry Smith     if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1908416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
190917ab2063SBarry Smith   }
191069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
191169957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1912bef8e0ddSBarry Smith 
1913bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
1914bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
1915bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
19163a40ed3dSBarry Smith   PetscFunctionReturn(0);
191717ab2063SBarry Smith }
191817ab2063SBarry Smith 
19195615d1e5SSatish Balay #undef __FUNC__
19205615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
19213d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
192217ab2063SBarry Smith {
1923416022c9SBarry Smith   Mat        C;
1924416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1925bef8e0ddSBarry Smith   int        i,len, m = a->m,shift = a->indexshift,ierr;
192617ab2063SBarry Smith 
19273a40ed3dSBarry Smith   PetscFunctionBegin;
19284043dd9cSLois Curfman McInnes   *B = 0;
1929f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView);
1930416022c9SBarry Smith   PLogObjectCreate(C);
19310452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1932f830108cSBarry Smith   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1933e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqAIJ;
1934e1311b90SBarry Smith   C->ops->view       = MatView_SeqAIJ;
1935416022c9SBarry Smith   C->factor     = A->factor;
1936416022c9SBarry Smith   c->row        = 0;
1937416022c9SBarry Smith   c->col        = 0;
193882bf6240SBarry Smith   c->icol       = 0;
1939416022c9SBarry Smith   c->indexshift = shift;
1940c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
194117ab2063SBarry Smith 
194244cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
194344cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
194444cd7ae7SLois Curfman McInnes   C->M          = a->m;
194544cd7ae7SLois Curfman McInnes   C->N          = a->n;
194617ab2063SBarry Smith 
19470452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
19480452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
194917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1950416022c9SBarry Smith     c->imax[i] = a->imax[i];
1951416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
195217ab2063SBarry Smith   }
195317ab2063SBarry Smith 
195417ab2063SBarry Smith   /* allocate the matrix space */
1955416022c9SBarry Smith   c->singlemalloc = 1;
1956416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
19570452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1958416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1959416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1960416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
196117ab2063SBarry Smith   if (m > 0) {
1962416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
196308480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1964416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
196517ab2063SBarry Smith     }
196608480c60SBarry Smith   }
196717ab2063SBarry Smith 
1968f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1969416022c9SBarry Smith   c->sorted      = a->sorted;
1970416022c9SBarry Smith   c->roworiented = a->roworiented;
1971416022c9SBarry Smith   c->nonew       = a->nonew;
19727a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
197317ab2063SBarry Smith 
1974416022c9SBarry Smith   if (a->diag) {
19750452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1976416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
197717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1978416022c9SBarry Smith       c->diag[i] = a->diag[i];
197917ab2063SBarry Smith     }
19803a40ed3dSBarry Smith   } else c->diag          = 0;
19816c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
19826c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1983754ec7b1SSatish Balay   if (a->inode.size){
1984daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1985754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1986daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1987754ec7b1SSatish Balay   } else {
1988754ec7b1SSatish Balay     c->inode.size       = 0;
1989754ec7b1SSatish Balay     c->inode.node_count = 0;
1990754ec7b1SSatish Balay   }
1991416022c9SBarry Smith   c->nz                 = a->nz;
1992416022c9SBarry Smith   c->maxnz              = a->maxnz;
1993416022c9SBarry Smith   c->solve_work         = 0;
199476dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1995754ec7b1SSatish Balay 
1996416022c9SBarry Smith   *B = C;
1997bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqAIJSetColumnIndices_C",
1998bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
1999bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
20003a40ed3dSBarry Smith   PetscFunctionReturn(0);
200117ab2063SBarry Smith }
200217ab2063SBarry Smith 
20035615d1e5SSatish Balay #undef __FUNC__
20045615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
200519bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
200617ab2063SBarry Smith {
2007416022c9SBarry Smith   Mat_SeqAIJ   *a;
2008416022c9SBarry Smith   Mat          B;
200917699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
2010bcd2baecSBarry Smith   MPI_Comm     comm;
201117ab2063SBarry Smith 
20123a40ed3dSBarry Smith   PetscFunctionBegin;
201319bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
201417699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
2015a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor");
201690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
20170752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
2018a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file");
201917ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
202017ab2063SBarry Smith 
2021d64ed03dSBarry Smith   if (nz < 0) {
2022a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ");
2023d64ed03dSBarry Smith   }
2024d64ed03dSBarry Smith 
202517ab2063SBarry Smith   /* read in row lengths */
20260452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
20270752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
202817ab2063SBarry Smith 
202917ab2063SBarry Smith   /* create our matrix */
2030416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
2031416022c9SBarry Smith   B = *A;
2032416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
2033416022c9SBarry Smith   shift = a->indexshift;
203417ab2063SBarry Smith 
203517ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
20360752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr);
203717ab2063SBarry Smith   if (shift) {
203817ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
2039416022c9SBarry Smith       a->j[i] += 1;
204017ab2063SBarry Smith     }
204117ab2063SBarry Smith   }
204217ab2063SBarry Smith 
204317ab2063SBarry Smith   /* read in nonzero values */
20440752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr);
204517ab2063SBarry Smith 
204617ab2063SBarry Smith   /* set matrix "i" values */
2047416022c9SBarry Smith   a->i[0] = -shift;
204817ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
2049416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2050416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
205117ab2063SBarry Smith   }
20520452661fSBarry Smith   PetscFree(rowlengths);
205317ab2063SBarry Smith 
20546d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20556d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20563a40ed3dSBarry Smith   PetscFunctionReturn(0);
205717ab2063SBarry Smith }
205817ab2063SBarry Smith 
20595615d1e5SSatish Balay #undef __FUNC__
20605615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
20618f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
20627264ac53SSatish Balay {
20637264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
20647264ac53SSatish Balay 
20653a40ed3dSBarry Smith   PetscFunctionBegin;
2066a8c6a408SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
20677264ac53SSatish Balay 
20687264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
20697264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
2070bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
20713a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2072bcd2baecSBarry Smith   }
20737264ac53SSatish Balay 
20747264ac53SSatish Balay   /* if the a->i are the same */
20758108c231SLois Curfman McInnes   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
20763a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
20777264ac53SSatish Balay   }
20787264ac53SSatish Balay 
20797264ac53SSatish Balay   /* if a->j are the same */
2080bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
20813a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2082bcd2baecSBarry Smith   }
2083bcd2baecSBarry Smith 
2084bcd2baecSBarry Smith   /* if a->a are the same */
208519bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
20863a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
20877264ac53SSatish Balay   }
208877c4ece6SBarry Smith   *flg = PETSC_TRUE;
20893a40ed3dSBarry Smith   PetscFunctionReturn(0);
20907264ac53SSatish Balay 
20917264ac53SSatish Balay }
2092