xref: /petsc/src/mat/impls/aij/seq/aij.c (revision cb5b572f1baa3b8b5e068cd6368cddce5382cbb3)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*cb5b572fSBarry Smith static char vcid[] = "$Id: aij.c,v 1.281 1998/08/20 15:37:03 bsmith Exp bsmith $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5d5d45c9bSBarry Smith /*
63369ce9aSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7d5d45c9bSBarry Smith   matrix storage format.
8d5d45c9bSBarry Smith */
93369ce9aSBarry Smith 
103369ce9aSBarry Smith #include "pinclude/pviewer.h"
113369ce9aSBarry Smith #include "sys.h"
1270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
13f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
14f5eb4b81SSatish Balay #include "src/inline/spops.h"
158d195f9aSBarry Smith #include "src/inline/dot.h"
16eeb667a2SSatish Balay #include "bitarray.h"
1717ab2063SBarry Smith 
18a2ce50c7SBarry Smith /*
19a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
20a2ce50c7SBarry Smith */
215615d1e5SSatish Balay #undef __FUNC__
225615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ"
23a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
24a2ce50c7SBarry Smith {
25a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
26a2ce50c7SBarry Smith   int        ierr = 1;
27a2ce50c7SBarry Smith 
283a40ed3dSBarry Smith   PetscFunctionBegin;
29e3372554SBarry Smith   SETERRQ(ierr,0,"Not implemented");
30d64ed03dSBarry Smith #if !defined(USE_PETSC_DEBUG)
31d64ed03dSBarry Smith   PetscFunctionReturn(0);
32d64ed03dSBarry Smith #endif
33a2ce50c7SBarry Smith }
34a2ce50c7SBarry Smith 
35bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3617ab2063SBarry Smith 
375615d1e5SSatish Balay #undef __FUNC__
385615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ"
398f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
406945ee14SBarry Smith                            PetscTruth *done)
4117ab2063SBarry Smith {
42416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
436945ee14SBarry Smith   int        ierr,i,ishift;
4417ab2063SBarry Smith 
453a40ed3dSBarry Smith   PetscFunctionBegin;
4631625ec5SSatish Balay   *m     = A->m;
473a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
486945ee14SBarry Smith   ishift = a->indexshift;
496945ee14SBarry Smith   if (symmetric) {
5031625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
516945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
5231625ec5SSatish Balay     int nz = a->i[a->m];
533b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
5431625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
553b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
563b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
5731625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
586945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
5931625ec5SSatish Balay     int nz = a->i[a->m] + 1;
603b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
6131625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
623b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
633b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
6431625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
656945ee14SBarry Smith   } else {
666945ee14SBarry Smith     *ia = a->i; *ja = a->j;
67a2ce50c7SBarry Smith   }
68a2ce50c7SBarry Smith 
693a40ed3dSBarry Smith   PetscFunctionReturn(0);
70a2744918SBarry Smith }
71a2744918SBarry Smith 
725615d1e5SSatish Balay #undef __FUNC__
735615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ"
748f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
756945ee14SBarry Smith                                PetscTruth *done)
766945ee14SBarry Smith {
776945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
783b2fbd54SBarry Smith   int        ishift = a->indexshift;
796945ee14SBarry Smith 
803a40ed3dSBarry Smith   PetscFunctionBegin;
813a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
823b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
836945ee14SBarry Smith     PetscFree(*ia);
846945ee14SBarry Smith     PetscFree(*ja);
85bcd2baecSBarry Smith   }
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
8717ab2063SBarry Smith }
8817ab2063SBarry Smith 
895615d1e5SSatish Balay #undef __FUNC__
905615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ"
9143a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
923b2fbd54SBarry Smith                            PetscTruth *done)
933b2fbd54SBarry Smith {
943b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
95a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
96a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
973b2fbd54SBarry Smith 
983a40ed3dSBarry Smith   PetscFunctionBegin;
993b2fbd54SBarry Smith   *nn     = A->n;
1003a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1013b2fbd54SBarry Smith   if (symmetric) {
102179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
1033b2fbd54SBarry Smith   } else {
10461d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
1053b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1063b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
107a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
1083b2fbd54SBarry Smith     jj = a->j;
1093b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
1103b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
1113b2fbd54SBarry Smith     }
1123b2fbd54SBarry Smith     cia[0] = oshift;
1133b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
1143b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
1153b2fbd54SBarry Smith     }
1163b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1173b2fbd54SBarry Smith     jj = a->j;
118a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
119a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
120a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1213b2fbd54SBarry Smith         col = *jj++ + ishift;
1223b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1233b2fbd54SBarry Smith       }
1243b2fbd54SBarry Smith     }
1253b2fbd54SBarry Smith     PetscFree(collengths);
1263b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1273b2fbd54SBarry Smith   }
1283b2fbd54SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1303b2fbd54SBarry Smith }
1313b2fbd54SBarry Smith 
1325615d1e5SSatish Balay #undef __FUNC__
1335615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ"
13443a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1353b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1363b2fbd54SBarry Smith {
1373a40ed3dSBarry Smith   PetscFunctionBegin;
1383a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1393b2fbd54SBarry Smith 
1403b2fbd54SBarry Smith   PetscFree(*ia);
1413b2fbd54SBarry Smith   PetscFree(*ja);
1423b2fbd54SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1443b2fbd54SBarry Smith }
1453b2fbd54SBarry Smith 
146227d817aSBarry Smith #define CHUNKSIZE   15
14717ab2063SBarry Smith 
1485615d1e5SSatish Balay #undef __FUNC__
1495615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ"
15061d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
15117ab2063SBarry Smith {
152416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
153416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1544b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
155d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
156416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
15717ab2063SBarry Smith 
1583a40ed3dSBarry Smith   PetscFunctionBegin;
15917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
160416022c9SBarry Smith     row  = im[k];
1613a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
162a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
163a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
1643b2fbd54SBarry Smith #endif
16517ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
16617ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
167416022c9SBarry Smith     low = 0;
16817ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
170a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
171a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
1723b2fbd54SBarry Smith #endif
1734b0e389bSBarry Smith       col = in[l] - shift;
1744b0e389bSBarry Smith       if (roworiented) {
1754b0e389bSBarry Smith         value = *v++;
176bef8e0ddSBarry Smith       } else {
1774b0e389bSBarry Smith         value = v[k + l*m];
1784b0e389bSBarry Smith       }
179416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
180416022c9SBarry Smith       while (high-low > 5) {
181416022c9SBarry Smith         t = (low+high)/2;
182416022c9SBarry Smith         if (rp[t] > col) high = t;
183416022c9SBarry Smith         else             low  = t;
18417ab2063SBarry Smith       }
185416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
18617ab2063SBarry Smith         if (rp[i] > col) break;
18717ab2063SBarry Smith         if (rp[i] == col) {
188416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
18917ab2063SBarry Smith           else                  ap[i] = value;
19017ab2063SBarry Smith           goto noinsert;
19117ab2063SBarry Smith         }
19217ab2063SBarry Smith       }
193c2653b3dSLois Curfman McInnes       if (nonew == 1) goto noinsert;
194a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
19517ab2063SBarry Smith       if (nrow >= rmax) {
19617ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
197416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
19817ab2063SBarry Smith         Scalar *new_a;
19917ab2063SBarry Smith 
200a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
20196854ed6SLois Curfman McInnes 
20217ab2063SBarry Smith         /* malloc new storage space */
203416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
2040452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
20517ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
20617ab2063SBarry Smith         new_i   = new_j + new_nz;
20717ab2063SBarry Smith 
20817ab2063SBarry Smith         /* copy over old data into new slots */
20917ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
210416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
211416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
212416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
213416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
21417ab2063SBarry Smith                                                            len*sizeof(int));
215416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
216416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
21717ab2063SBarry Smith                                                            len*sizeof(Scalar));
21817ab2063SBarry Smith         /* free up old matrix storage */
2190452661fSBarry Smith         PetscFree(a->a);
2200452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
221416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
222416022c9SBarry Smith         a->singlemalloc = 1;
22317ab2063SBarry Smith 
22417ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
225416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
226416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
227416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
228b810aeb4SBarry Smith         a->reallocs++;
22917ab2063SBarry Smith       }
230416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
231416022c9SBarry Smith       /* shift up all the later entries in this row */
232416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
23317ab2063SBarry Smith         rp[ii+1] = rp[ii];
23417ab2063SBarry Smith         ap[ii+1] = ap[ii];
23517ab2063SBarry Smith       }
23617ab2063SBarry Smith       rp[i] = col;
23717ab2063SBarry Smith       ap[i] = value;
23817ab2063SBarry Smith       noinsert:;
239416022c9SBarry Smith       low = i + 1;
24017ab2063SBarry Smith     }
24117ab2063SBarry Smith     ailen[row] = nrow;
24217ab2063SBarry Smith   }
2433a40ed3dSBarry Smith   PetscFunctionReturn(0);
24417ab2063SBarry Smith }
24517ab2063SBarry Smith 
2465615d1e5SSatish Balay #undef __FUNC__
2475615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ"
2488f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2497eb43aa7SLois Curfman McInnes {
2507eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
251b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2527eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2537eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2547eb43aa7SLois Curfman McInnes 
2553a40ed3dSBarry Smith   PetscFunctionBegin;
2567eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2577eb43aa7SLois Curfman McInnes     row  = im[k];
258a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
259a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
2607eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2617eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2627eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
263a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
264a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
2657eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2667eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2677eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2687eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2697eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2707eb43aa7SLois Curfman McInnes         else             low  = t;
2717eb43aa7SLois Curfman McInnes       }
2727eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2737eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2747eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
275b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2767eb43aa7SLois Curfman McInnes           goto finished;
2777eb43aa7SLois Curfman McInnes         }
2787eb43aa7SLois Curfman McInnes       }
279b49de8d1SLois Curfman McInnes       *v++ = zero;
2807eb43aa7SLois Curfman McInnes       finished:;
2817eb43aa7SLois Curfman McInnes     }
2827eb43aa7SLois Curfman McInnes   }
2833a40ed3dSBarry Smith   PetscFunctionReturn(0);
2847eb43aa7SLois Curfman McInnes }
2857eb43aa7SLois Curfman McInnes 
28617ab2063SBarry Smith 
2875615d1e5SSatish Balay #undef __FUNC__
2885615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary"
289480ef9eaSBarry Smith int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
29017ab2063SBarry Smith {
291416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
292416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
29317ab2063SBarry Smith 
2943a40ed3dSBarry Smith   PetscFunctionBegin;
29590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2960452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
297416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
298416022c9SBarry Smith   col_lens[1] = a->m;
299416022c9SBarry Smith   col_lens[2] = a->n;
300416022c9SBarry Smith   col_lens[3] = a->nz;
301416022c9SBarry Smith 
302416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
303416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
304416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
30517ab2063SBarry Smith   }
3060752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3070452661fSBarry Smith   PetscFree(col_lens);
308416022c9SBarry Smith 
309416022c9SBarry Smith   /* store column indices (zero start index) */
310416022c9SBarry Smith   if (a->indexshift) {
311416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
31217ab2063SBarry Smith   }
3130752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0); CHKERRQ(ierr);
314416022c9SBarry Smith   if (a->indexshift) {
315416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
31617ab2063SBarry Smith   }
317416022c9SBarry Smith 
318416022c9SBarry Smith   /* store nonzero values */
3190752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
32117ab2063SBarry Smith }
322416022c9SBarry Smith 
3235615d1e5SSatish Balay #undef __FUNC__
3245615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII"
325480ef9eaSBarry Smith int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
326416022c9SBarry Smith {
327416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
328496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
32917ab2063SBarry Smith   FILE        *fd;
33017ab2063SBarry Smith   char        *outputname;
33117ab2063SBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
33390ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
334416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
33590ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
336a93ec695SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
3373a40ed3dSBarry Smith     PetscFunctionReturn(0);
3383a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
339496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
340496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
341496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
34295e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
34395e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
3443a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
345d00d2cf4SBarry Smith     int nofinalvalue = 0;
346d00d2cf4SBarry Smith     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) {
347d00d2cf4SBarry Smith       nofinalvalue = 1;
348d00d2cf4SBarry Smith     }
349416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3504e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
351d00d2cf4SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue);
35217ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
35317ab2063SBarry Smith 
35417ab2063SBarry Smith     for (i=0; i<m; i++) {
355416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
3563a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
357e20fef11SSatish Balay         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
35817ab2063SBarry Smith #else
3597a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
36017ab2063SBarry Smith #endif
36117ab2063SBarry Smith       }
36217ab2063SBarry Smith     }
363d00d2cf4SBarry Smith     if (nofinalvalue) {
364d00d2cf4SBarry Smith       fprintf(fd,"%d %d  %18.16e\n", m, a->n, 0.0);
365d00d2cf4SBarry Smith     }
36617ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
3673a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
36844cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
36944cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
37044cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
3713a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
372e20fef11SSatish Balay         if (PetscImaginary(a->a[j]) > 0.0 && PetscReal(a->a[j]) != 0.0)
373e20fef11SSatish Balay           fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
374e20fef11SSatish Balay         else if (PetscImaginary(a->a[j]) < 0.0 && PetscReal(a->a[j]) != 0.0)
375e20fef11SSatish Balay           fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j]));
376e20fef11SSatish Balay         else if (PetscReal(a->a[j]) != 0.0)
377e20fef11SSatish Balay           fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j]));
37844cd7ae7SLois Curfman McInnes #else
37944cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
38044cd7ae7SLois Curfman McInnes #endif
38144cd7ae7SLois Curfman McInnes       }
38244cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
38344cd7ae7SLois Curfman McInnes     }
38444cd7ae7SLois Curfman McInnes   }
385496be53dSLois Curfman McInnes   else if (format == VIEWER_FORMAT_ASCII_SYMMODU) {
386496be53dSLois Curfman McInnes     int nzd=0, fshift=1, *sptr;
3872e44a96cSLois Curfman McInnes     sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr);
388496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
389496be53dSLois Curfman McInnes       sptr[i] = nzd+1;
390496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
391496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
3923a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
393e20fef11SSatish Balay           if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) nzd++;
394496be53dSLois Curfman McInnes #else
395496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) nzd++;
396496be53dSLois Curfman McInnes #endif
397496be53dSLois Curfman McInnes         }
398496be53dSLois Curfman McInnes       }
399496be53dSLois Curfman McInnes     }
4002e44a96cSLois Curfman McInnes     sptr[m] = nzd+1;
401496be53dSLois Curfman McInnes     fprintf(fd," %d %d\n\n",m,nzd);
4022e44a96cSLois Curfman McInnes     for ( i=0; i<m+1; i+=6 ) {
4032e44a96cSLois Curfman McInnes       if (i+4<m) fprintf(fd," %d %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);
4042e44a96cSLois Curfman McInnes       else if (i+3<m) fprintf(fd," %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);
4052e44a96cSLois Curfman McInnes       else if (i+2<m) fprintf(fd," %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);
4062e44a96cSLois Curfman McInnes       else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]);
4072e44a96cSLois Curfman McInnes       else if (i<m)   fprintf(fd," %d %d\n",sptr[i],sptr[i+1]);
4087272d637SLois Curfman McInnes       else            fprintf(fd," %d\n",sptr[i]);
409496be53dSLois Curfman McInnes     }
410496be53dSLois Curfman McInnes     fprintf(fd,"\n");
411496be53dSLois Curfman McInnes     PetscFree(sptr);
412496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
413496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
414496be53dSLois Curfman McInnes         if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift);
415496be53dSLois Curfman McInnes       }
416496be53dSLois Curfman McInnes       fprintf(fd,"\n");
417496be53dSLois Curfman McInnes     }
418496be53dSLois Curfman McInnes     fprintf(fd,"\n");
419496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
420496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
421496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
4223a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
423e20fef11SSatish Balay           if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0)
424e20fef11SSatish Balay             fprintf(fd," %18.16e %18.16e ",PetscReal(a->a[j]),PetscImaginary(a->a[j]));
425496be53dSLois Curfman McInnes #else
426496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]);
427496be53dSLois Curfman McInnes #endif
428496be53dSLois Curfman McInnes         }
429496be53dSLois Curfman McInnes       }
430496be53dSLois Curfman McInnes       fprintf(fd,"\n");
431496be53dSLois Curfman McInnes     }
4323a40ed3dSBarry Smith   } else {
43317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
43417ab2063SBarry Smith       fprintf(fd,"row %d:",i);
435416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
4363a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
437e20fef11SSatish Balay         if (PetscImaginary(a->a[j]) > 0.0) {
438e20fef11SSatish Balay           fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j]));
439e20fef11SSatish Balay         } else if (PetscImaginary(a->a[j]) < 0.0) {
440e20fef11SSatish Balay           fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j]));
4413a40ed3dSBarry Smith         } else {
442e20fef11SSatish Balay           fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j]));
44317ab2063SBarry Smith         }
44417ab2063SBarry Smith #else
445416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
44617ab2063SBarry Smith #endif
44717ab2063SBarry Smith       }
44817ab2063SBarry Smith       fprintf(fd,"\n");
44917ab2063SBarry Smith     }
45017ab2063SBarry Smith   }
45117ab2063SBarry Smith   fflush(fd);
4523a40ed3dSBarry Smith   PetscFunctionReturn(0);
453416022c9SBarry Smith }
454416022c9SBarry Smith 
4555615d1e5SSatish Balay #undef __FUNC__
456480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw_Zoom"
457480ef9eaSBarry Smith int MatView_SeqAIJ_Draw_Zoom(Draw draw,void *Aa)
458416022c9SBarry Smith {
459480ef9eaSBarry Smith   Mat         A = (Mat) Aa;
460416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
461480ef9eaSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,color;
4620513a670SBarry Smith   int         format;
463480ef9eaSBarry Smith   double      xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
464480ef9eaSBarry Smith   Viewer      viewer;
465cddf8d76SBarry Smith 
4663a40ed3dSBarry Smith   PetscFunctionBegin;
467480ef9eaSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr);
4680513a670SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
46919bcc07fSBarry Smith 
470480ef9eaSBarry Smith   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr);
471416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
4720513a670SBarry Smith 
4730513a670SBarry Smith   if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4740513a670SBarry Smith     /* Blue for negative, Cyan for zero and  Red for positive */
475cddf8d76SBarry Smith     color = DRAW_BLUE;
476416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
477cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
478416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
479cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
4803a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
481e20fef11SSatish Balay         if (PetscReal(a->a[j]) >=  0.) continue;
482cddf8d76SBarry Smith #else
483cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
484cddf8d76SBarry Smith #endif
485480ef9eaSBarry Smith         ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
486cddf8d76SBarry Smith       }
487cddf8d76SBarry Smith     }
488cddf8d76SBarry Smith     color = DRAW_CYAN;
489cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
490cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
491cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
492cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
493cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
494480ef9eaSBarry Smith         ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
495cddf8d76SBarry Smith       }
496cddf8d76SBarry Smith     }
497cddf8d76SBarry Smith     color = DRAW_RED;
498cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
499cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
500cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
501cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5023a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
503e20fef11SSatish Balay         if (PetscReal(a->a[j]) <=  0.) continue;
504cddf8d76SBarry Smith #else
505cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
506cddf8d76SBarry Smith #endif
507480ef9eaSBarry Smith         ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
508416022c9SBarry Smith       }
509416022c9SBarry Smith     }
5100513a670SBarry Smith   } else {
5110513a670SBarry Smith     /* use contour shading to indicate magnitude of values */
5120513a670SBarry Smith     /* first determine max of all nonzero values */
5130513a670SBarry Smith     int    nz = a->nz,count;
5140513a670SBarry Smith     Draw   popup;
515480ef9eaSBarry Smith     double scale;
5160513a670SBarry Smith 
5170513a670SBarry Smith     for ( i=0; i<nz; i++ ) {
5180513a670SBarry Smith       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
5190513a670SBarry Smith     }
520480ef9eaSBarry Smith     scale = (245.0 - DRAW_BASIC_COLORS)/maxv;
521522c5e43SBarry Smith     ierr  = DrawGetPopup(draw,&popup); CHKERRQ(ierr);
5220513a670SBarry Smith     ierr  = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr);
5230513a670SBarry Smith     count = 0;
5240513a670SBarry Smith     for ( i=0; i<m; i++ ) {
5250513a670SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
5260513a670SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5270513a670SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5286d6b6c30SSatish Balay         color = DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(a->a[count]));
529480ef9eaSBarry Smith         ierr  = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
5300513a670SBarry Smith         count++;
5310513a670SBarry Smith       }
5320513a670SBarry Smith     }
5330513a670SBarry Smith   }
534480ef9eaSBarry Smith   PetscFunctionReturn(0);
535480ef9eaSBarry Smith }
536cddf8d76SBarry Smith 
537480ef9eaSBarry Smith #undef __FUNC__
538480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw"
539480ef9eaSBarry Smith int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
540480ef9eaSBarry Smith {
541480ef9eaSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data;
542480ef9eaSBarry Smith   int        ierr;
543480ef9eaSBarry Smith   Draw       draw;
544480ef9eaSBarry Smith   double     xr,yr,xl,yl,h,w;
545480ef9eaSBarry Smith 
546480ef9eaSBarry Smith   PetscTruth isnull;
547480ef9eaSBarry Smith 
548480ef9eaSBarry Smith   PetscFunctionBegin;
549480ef9eaSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
550480ef9eaSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr);
551480ef9eaSBarry Smith   if (isnull) PetscFunctionReturn(0);
552480ef9eaSBarry Smith 
553480ef9eaSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
554480ef9eaSBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
555480ef9eaSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
556cddf8d76SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
557480ef9eaSBarry Smith   ierr = DrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A); CHKERRQ(ierr);
558480ef9eaSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5593a40ed3dSBarry Smith   PetscFunctionReturn(0);
560416022c9SBarry Smith }
561416022c9SBarry Smith 
5625615d1e5SSatish Balay #undef __FUNC__
563d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ"
564e1311b90SBarry Smith int MatView_SeqAIJ(Mat A,Viewer viewer)
565416022c9SBarry Smith {
566416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
567bcd2baecSBarry Smith   ViewerType  vtype;
568bcd2baecSBarry Smith   int         ierr;
569416022c9SBarry Smith 
5703a40ed3dSBarry Smith   PetscFunctionBegin;
571bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
572bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
5733a40ed3dSBarry Smith     ierr = ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);  CHKERRQ(ierr);
5745cd90555SBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
5753a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_ASCII(A,viewer); CHKERRQ(ierr);
5765cd90555SBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
5773a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_Binary(A,viewer); CHKERRQ(ierr);
5785cd90555SBarry Smith   } else if (vtype == DRAW_VIEWER) {
5793a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_Draw(A,viewer); CHKERRQ(ierr);
5805cd90555SBarry Smith   } else {
5815cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
58217ab2063SBarry Smith   }
5833a40ed3dSBarry Smith   PetscFunctionReturn(0);
58417ab2063SBarry Smith }
58519bcc07fSBarry Smith 
586c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
5875615d1e5SSatish Balay #undef __FUNC__
5885615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ"
5898f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
59017ab2063SBarry Smith {
591416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
59241c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
59343ee02c3SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0;
594416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
59517ab2063SBarry Smith 
5963a40ed3dSBarry Smith   PetscFunctionBegin;
5973a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
59817ab2063SBarry Smith 
59943ee02c3SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
60017ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
601416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
60217ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
60394a9d846SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
60417ab2063SBarry Smith     if (fshift) {
605416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
60617ab2063SBarry Smith       N = ailen[i];
60717ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
60817ab2063SBarry Smith         ip[j-fshift] = ip[j];
60917ab2063SBarry Smith         ap[j-fshift] = ap[j];
61017ab2063SBarry Smith       }
61117ab2063SBarry Smith     }
61217ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
61317ab2063SBarry Smith   }
61417ab2063SBarry Smith   if (m) {
61517ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
61617ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
61717ab2063SBarry Smith   }
61817ab2063SBarry Smith   /* reset ilen and imax for each row */
61917ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
62017ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
62117ab2063SBarry Smith   }
622416022c9SBarry Smith   a->nz = ai[m] + shift;
62317ab2063SBarry Smith 
62417ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
625416022c9SBarry Smith   if (fshift && a->diag) {
6260452661fSBarry Smith     PetscFree(a->diag);
627416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
628416022c9SBarry Smith     a->diag = 0;
62917ab2063SBarry Smith   }
6304e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
6314e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
6324e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
633b810aeb4SBarry Smith            a->reallocs);
63494a9d846SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax);
635dd5f02e7SSatish Balay   a->reallocs          = 0;
6364e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
6374e220ebcSLois Curfman McInnes 
63876dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
63941c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
6403a40ed3dSBarry Smith   PetscFunctionReturn(0);
64117ab2063SBarry Smith }
64217ab2063SBarry Smith 
6435615d1e5SSatish Balay #undef __FUNC__
6445615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ"
6458f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A)
64617ab2063SBarry Smith {
647416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
6483a40ed3dSBarry Smith 
6493a40ed3dSBarry Smith   PetscFunctionBegin;
650cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
6513a40ed3dSBarry Smith   PetscFunctionReturn(0);
65217ab2063SBarry Smith }
653416022c9SBarry Smith 
6545615d1e5SSatish Balay #undef __FUNC__
6555615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ"
656e1311b90SBarry Smith int MatDestroy_SeqAIJ(Mat A)
65717ab2063SBarry Smith {
658416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
65982bf6240SBarry Smith   int        ierr;
660d5d45c9bSBarry Smith 
6613a40ed3dSBarry Smith   PetscFunctionBegin;
66294d884c6SBarry Smith   if (--A->refct > 0) PetscFunctionReturn(0);
66394d884c6SBarry Smith 
66494d884c6SBarry Smith   if (A->mapping) {
66594d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr);
66694d884c6SBarry Smith   }
66794d884c6SBarry Smith   if (A->bmapping) {
66894d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr);
66994d884c6SBarry Smith   }
67061b13de0SBarry Smith   if (A->rmap) {
67161b13de0SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
67261b13de0SBarry Smith   }
67361b13de0SBarry Smith   if (A->cmap) {
67461b13de0SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
67561b13de0SBarry Smith   }
6763a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
677e1311b90SBarry Smith   PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
67817ab2063SBarry Smith #endif
6790452661fSBarry Smith   PetscFree(a->a);
6800452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
6810452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
6820452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
6830452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
6840452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
68576dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
68682bf6240SBarry Smith   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
6870452661fSBarry Smith   PetscFree(a);
688eed86810SBarry Smith 
689f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
690f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
6913a40ed3dSBarry Smith   PetscFunctionReturn(0);
69217ab2063SBarry Smith }
69317ab2063SBarry Smith 
6945615d1e5SSatish Balay #undef __FUNC__
6955615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
6968f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A)
69717ab2063SBarry Smith {
6983a40ed3dSBarry Smith   PetscFunctionBegin;
6993a40ed3dSBarry Smith   PetscFunctionReturn(0);
70017ab2063SBarry Smith }
70117ab2063SBarry Smith 
7025615d1e5SSatish Balay #undef __FUNC__
7035615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
7048f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op)
70517ab2063SBarry Smith {
706416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7073a40ed3dSBarry Smith 
7083a40ed3dSBarry Smith   PetscFunctionBegin;
7096d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
7106d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
7116d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
712219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
7136d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
714c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
71596854ed6SLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
7166d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
7176d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
718219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7196d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7206d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
72190f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
722b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES||
723b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
724981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n");
7253a40ed3dSBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS) {
7263a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7273a40ed3dSBarry Smith   } else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
7286d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
7296d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
7306d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
7316d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
7323a40ed3dSBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7333a40ed3dSBarry Smith   PetscFunctionReturn(0);
73417ab2063SBarry Smith }
73517ab2063SBarry Smith 
7365615d1e5SSatish Balay #undef __FUNC__
7375615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
7388f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
73917ab2063SBarry Smith {
740416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7413a40ed3dSBarry Smith   int        i,j, n,shift = a->indexshift,ierr;
74217ab2063SBarry Smith   Scalar     *x, zero = 0.0;
74317ab2063SBarry Smith 
7443a40ed3dSBarry Smith   PetscFunctionBegin;
7453a40ed3dSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
746e1311b90SBarry Smith   ierr = VecGetArray(v,&x);;CHKERRQ(ierr);
747e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr);
748a8c6a408SBarry Smith   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
749416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
750416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
751416022c9SBarry Smith       if (a->j[j]+shift == i) {
752416022c9SBarry Smith         x[i] = a->a[j];
75317ab2063SBarry Smith         break;
75417ab2063SBarry Smith       }
75517ab2063SBarry Smith     }
75617ab2063SBarry Smith   }
757e1311b90SBarry Smith   ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr);
7583a40ed3dSBarry Smith   PetscFunctionReturn(0);
75917ab2063SBarry Smith }
76017ab2063SBarry Smith 
76117ab2063SBarry Smith /* -------------------------------------------------------*/
76217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
76317ab2063SBarry Smith /* -------------------------------------------------------*/
7645615d1e5SSatish Balay #undef __FUNC__
7655615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
76644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
76717ab2063SBarry Smith {
768416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
76917ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
770e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx, shift = a->indexshift;
77117ab2063SBarry Smith 
7723a40ed3dSBarry Smith   PetscFunctionBegin;
773e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
774e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
775cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
77617ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
77717ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
778416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
779416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
780416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
78117ab2063SBarry Smith     alpha = x[i];
78217ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
78317ab2063SBarry Smith   }
784416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
785e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
786e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7873a40ed3dSBarry Smith   PetscFunctionReturn(0);
78817ab2063SBarry Smith }
789d5d45c9bSBarry Smith 
7905615d1e5SSatish Balay #undef __FUNC__
7915615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
79244cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
79317ab2063SBarry Smith {
794416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
79517ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
796e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx,shift = a->indexshift;
79717ab2063SBarry Smith 
7983a40ed3dSBarry Smith   PetscFunctionBegin;
799e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
800e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
80117ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
80217ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
80317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
804416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
805416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
806416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
80717ab2063SBarry Smith     alpha = x[i];
80817ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
80917ab2063SBarry Smith   }
81090f02eecSBarry Smith   PLogFlops(2*a->nz);
811e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
812e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8133a40ed3dSBarry Smith   PetscFunctionReturn(0);
81417ab2063SBarry Smith }
81517ab2063SBarry Smith 
8165615d1e5SSatish Balay #undef __FUNC__
8175615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
81844cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
81917ab2063SBarry Smith {
820416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
82117ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
822e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
8232e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTAIJ)
824e36a17ebSSatish Balay   int        n, i, jrow,j;
825e36a17ebSSatish Balay #endif
82617ab2063SBarry Smith 
827fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT)
828fee21e36SBarry Smith #pragma disjoint(*x,*y,*v)
829fee21e36SBarry Smith #endif
830fee21e36SBarry Smith 
8313a40ed3dSBarry Smith   PetscFunctionBegin;
832e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
833e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
83417ab2063SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
835416022c9SBarry Smith   idx  = a->j;
836d64ed03dSBarry Smith   v    = a->a;
837416022c9SBarry Smith   ii   = a->i;
838acc4d009SSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTAIJ)
83929b5ca8aSSatish Balay   fortranmultaij_(&m,x,ii,idx+shift,v+shift,y);
8408d195f9aSBarry Smith #else
841d64ed03dSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
842d64ed03dSBarry Smith   idx  += shift;
84317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8449ea0dfa2SSatish Balay     jrow = ii[i];
8459ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
84617ab2063SBarry Smith     sum  = 0.0;
8479ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8489ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8499ea0dfa2SSatish Balay      }
85017ab2063SBarry Smith     y[i] = sum;
85117ab2063SBarry Smith   }
8528d195f9aSBarry Smith #endif
853416022c9SBarry Smith   PLogFlops(2*a->nz - m);
854e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
855e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
8563a40ed3dSBarry Smith   PetscFunctionReturn(0);
85717ab2063SBarry Smith }
85817ab2063SBarry Smith 
8595615d1e5SSatish Balay #undef __FUNC__
8605615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
86144cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
86217ab2063SBarry Smith {
863416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
86417ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
865e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
8662e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTADDAIJ)
867e36a17ebSSatish Balay   int        n,i,jrow,j;
868e36a17ebSSatish Balay #endif
8699ea0dfa2SSatish Balay 
8703a40ed3dSBarry Smith   PetscFunctionBegin;
871e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
872e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
873e1311b90SBarry Smith   ierr = VecGetArray(zz,&z); CHKERRQ(ierr);
87417ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
875cddf8d76SBarry Smith   idx  = a->j;
876d64ed03dSBarry Smith   v    = a->a;
877cddf8d76SBarry Smith   ii   = a->i;
8782e40c62eSSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTADDAIJ)
87929b5ca8aSSatish Balay   fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
88002ab625aSSatish Balay #else
881d64ed03dSBarry Smith   v   += shift; /* shift for Fortran start by 1 indexing */
882d64ed03dSBarry Smith   idx += shift;
88317ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8849ea0dfa2SSatish Balay     jrow = ii[i];
8859ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
88617ab2063SBarry Smith     sum  = y[i];
8879ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8889ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8899ea0dfa2SSatish Balay      }
89017ab2063SBarry Smith     z[i] = sum;
89117ab2063SBarry Smith   }
89202ab625aSSatish Balay #endif
893416022c9SBarry Smith   PLogFlops(2*a->nz);
894e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
895e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
896e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr);
8973a40ed3dSBarry Smith   PetscFunctionReturn(0);
89817ab2063SBarry Smith }
89917ab2063SBarry Smith 
90017ab2063SBarry Smith /*
90117ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
90217ab2063SBarry Smith */
90317ab2063SBarry Smith 
9045615d1e5SSatish Balay #undef __FUNC__
9055615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
906416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
90717ab2063SBarry Smith {
908416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
909416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
91017ab2063SBarry Smith 
9113a40ed3dSBarry Smith   PetscFunctionBegin;
9120452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
913416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
914416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
91535b0346bSBarry Smith     diag[i] = a->i[i+1];
916416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
917416022c9SBarry Smith       if (a->j[j]+shift == i) {
91817ab2063SBarry Smith         diag[i] = j - shift;
91917ab2063SBarry Smith         break;
92017ab2063SBarry Smith       }
92117ab2063SBarry Smith     }
92217ab2063SBarry Smith   }
923416022c9SBarry Smith   a->diag = diag;
9243a40ed3dSBarry Smith   PetscFunctionReturn(0);
92517ab2063SBarry Smith }
92617ab2063SBarry Smith 
9275615d1e5SSatish Balay #undef __FUNC__
9285615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
92944cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
93017ab2063SBarry Smith                            double fshift,int its,Vec xx)
93117ab2063SBarry Smith {
932416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
933416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
934d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
93517ab2063SBarry Smith 
9363a40ed3dSBarry Smith   PetscFunctionBegin;
937e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
938e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
939d00d2cf4SBarry Smith   if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);}
940416022c9SBarry Smith   diag = a->diag;
941416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
94217ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
94317ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
94417ab2063SBarry Smith     bs = b + shift;
94517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
946416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
947416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
948488ecbafSBarry Smith 	PLogFlops(2*n-1);
949416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
950416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
95117ab2063SBarry Smith         sum  = b[i]*d/omega;
95217ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
95317ab2063SBarry Smith         x[i] = sum;
95417ab2063SBarry Smith     }
955*cb5b572fSBarry Smith     ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
956*cb5b572fSBarry Smith     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
9573a40ed3dSBarry Smith     PetscFunctionReturn(0);
95817ab2063SBarry Smith   }
95917ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
960a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done");
9613a40ed3dSBarry Smith   } else if (flag & SOR_EISENSTAT) {
96217ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
96317ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
96417ab2063SBarry Smith 
96517ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
96617ab2063SBarry Smith 
96717ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
96817ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
96917ab2063SBarry Smith     is the relaxation factor.
97017ab2063SBarry Smith     */
9710452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
97217ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
97317ab2063SBarry Smith 
97417ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
97517ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
976416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
977416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
978488ecbafSBarry Smith       PLogFlops(2*n-1);
979416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
980416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
98117ab2063SBarry Smith       sum  = b[i];
98217ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
98317ab2063SBarry Smith       x[i] = omega*(sum/d);
98417ab2063SBarry Smith     }
98517ab2063SBarry Smith 
98617ab2063SBarry Smith     /*  t = b - (2*E - D)x */
987416022c9SBarry Smith     v = a->a;
98817ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
989488ecbafSBarry Smith     PLogFlops(2*m);
990488ecbafSBarry Smith 
99117ab2063SBarry Smith 
99217ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
993416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
994416022c9SBarry Smith     diag = a->diag;
99517ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
996416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
997416022c9SBarry Smith       n    = diag[i] - a->i[i];
998488ecbafSBarry Smith       PLogFlops(2*n-1);
999416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
1000416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
100117ab2063SBarry Smith       sum  = t[i];
100217ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
100317ab2063SBarry Smith       t[i] = omega*(sum/d);
100417ab2063SBarry Smith     }
100517ab2063SBarry Smith 
100617ab2063SBarry Smith     /*  x = x + t */
100717ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
10080452661fSBarry Smith     PetscFree(t);
1009*cb5b572fSBarry Smith     ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
1010*cb5b572fSBarry Smith     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10113a40ed3dSBarry Smith     PetscFunctionReturn(0);
101217ab2063SBarry Smith   }
101317ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
101417ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
101517ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1016416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1017416022c9SBarry Smith         n    = diag[i] - a->i[i];
1018488ecbafSBarry Smith 	PLogFlops(2*n-1);
1019416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1020416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
102117ab2063SBarry Smith         sum  = b[i];
102217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
102317ab2063SBarry Smith         x[i] = omega*(sum/d);
102417ab2063SBarry Smith       }
102517ab2063SBarry Smith       xb = x;
10263a40ed3dSBarry Smith     } else xb = b;
102717ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
102817ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
102917ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1030416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
103117ab2063SBarry Smith       }
1032488ecbafSBarry Smith       PLogFlops(m);
103317ab2063SBarry Smith     }
103417ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
103517ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1036416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1037416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
1038488ecbafSBarry Smith 	PLogFlops(2*n-1);
1039416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1040416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
104117ab2063SBarry Smith         sum  = xb[i];
104217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
104317ab2063SBarry Smith         x[i] = omega*(sum/d);
104417ab2063SBarry Smith       }
104517ab2063SBarry Smith     }
104617ab2063SBarry Smith     its--;
104717ab2063SBarry Smith   }
104817ab2063SBarry Smith   while (its--) {
104917ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
105017ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1051416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1052416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1053488ecbafSBarry Smith 	PLogFlops(2*n-1);
1054416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1055416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
105617ab2063SBarry Smith         sum  = b[i];
105717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10587e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
105917ab2063SBarry Smith       }
106017ab2063SBarry Smith     }
106117ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
106217ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1063416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1064416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1065488ecbafSBarry Smith 	PLogFlops(2*n-1);
1066416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1067416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
106817ab2063SBarry Smith         sum  = b[i];
106917ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10707e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
107117ab2063SBarry Smith       }
107217ab2063SBarry Smith     }
107317ab2063SBarry Smith   }
1074*cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
1075*cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10763a40ed3dSBarry Smith   PetscFunctionReturn(0);
107717ab2063SBarry Smith }
107817ab2063SBarry Smith 
10795615d1e5SSatish Balay #undef __FUNC__
10805615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10818f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
108217ab2063SBarry Smith {
1083416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
10844e220ebcSLois Curfman McInnes 
10853a40ed3dSBarry Smith   PetscFunctionBegin;
10864e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
10874e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
10884e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
10894e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
10904e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
10914e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
10924e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
10934e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
10944e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
10954e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
10964e220ebcSLois Curfman McInnes   info->memory         = A->mem;
10974e220ebcSLois Curfman McInnes   if (A->factor) {
10984e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
10994e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
11004e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
11014e220ebcSLois Curfman McInnes   } else {
11024e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
11034e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
11044e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
11054e220ebcSLois Curfman McInnes   }
11063a40ed3dSBarry Smith   PetscFunctionReturn(0);
110717ab2063SBarry Smith }
110817ab2063SBarry Smith 
110917ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
111017ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
111117ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
111217ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
111317ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
111417ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
111517ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
111617ab2063SBarry Smith 
11175615d1e5SSatish Balay #undef __FUNC__
11185615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
11198f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
112017ab2063SBarry Smith {
1121416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1122416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
112317ab2063SBarry Smith 
11243a40ed3dSBarry Smith   PetscFunctionBegin;
112577c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
112617ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
112717ab2063SBarry Smith   if (diag) {
112817ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1129a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1130416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1131416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1132416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1133416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
11343a40ed3dSBarry Smith       } else {
1135d64ed03dSBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr);
113617ab2063SBarry Smith       }
113717ab2063SBarry Smith     }
11383a40ed3dSBarry Smith   } else {
113917ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1140a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1141416022c9SBarry Smith       a->ilen[rows[i]] = 0;
114217ab2063SBarry Smith     }
114317ab2063SBarry Smith   }
114417ab2063SBarry Smith   ISRestoreIndices(is,&rows);
114543a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11463a40ed3dSBarry Smith   PetscFunctionReturn(0);
114717ab2063SBarry Smith }
114817ab2063SBarry Smith 
11495615d1e5SSatish Balay #undef __FUNC__
11505615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11518f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
115217ab2063SBarry Smith {
1153416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11543a40ed3dSBarry Smith 
11553a40ed3dSBarry Smith   PetscFunctionBegin;
11560752156aSBarry Smith   if (m) *m = a->m;
11570752156aSBarry Smith   if (n) *n = a->n;
11583a40ed3dSBarry Smith   PetscFunctionReturn(0);
115917ab2063SBarry Smith }
116017ab2063SBarry Smith 
11615615d1e5SSatish Balay #undef __FUNC__
11625615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11638f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
116417ab2063SBarry Smith {
1165416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11663a40ed3dSBarry Smith 
11673a40ed3dSBarry Smith   PetscFunctionBegin;
1168416022c9SBarry Smith   *m = 0; *n = a->m;
11693a40ed3dSBarry Smith   PetscFunctionReturn(0);
117017ab2063SBarry Smith }
1171026e39d0SSatish Balay 
11725615d1e5SSatish Balay #undef __FUNC__
11735615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11744e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
117517ab2063SBarry Smith {
1176416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1177c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
117817ab2063SBarry Smith 
11793a40ed3dSBarry Smith   PetscFunctionBegin;
1180a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
118117ab2063SBarry Smith 
1182416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1183416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
118417ab2063SBarry Smith   if (idx) {
1185416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
11864e093b46SBarry Smith     if (*nz && shift) {
11870452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
118817ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
11894e093b46SBarry Smith     } else if (*nz) {
11904e093b46SBarry Smith       *idx = itmp;
119117ab2063SBarry Smith     }
119217ab2063SBarry Smith     else *idx = 0;
119317ab2063SBarry Smith   }
11943a40ed3dSBarry Smith   PetscFunctionReturn(0);
119517ab2063SBarry Smith }
119617ab2063SBarry Smith 
11975615d1e5SSatish Balay #undef __FUNC__
11985615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
11994e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
120017ab2063SBarry Smith {
12014e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
12023a40ed3dSBarry Smith 
12033a40ed3dSBarry Smith   PetscFunctionBegin;
12044e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
12053a40ed3dSBarry Smith   PetscFunctionReturn(0);
120617ab2063SBarry Smith }
120717ab2063SBarry Smith 
12085615d1e5SSatish Balay #undef __FUNC__
12095615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
12108f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
121117ab2063SBarry Smith {
1212416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1213416022c9SBarry Smith   Scalar     *v = a->a;
121417ab2063SBarry Smith   double     sum = 0.0;
1215416022c9SBarry Smith   int        i, j,shift = a->indexshift;
121617ab2063SBarry Smith 
12173a40ed3dSBarry Smith   PetscFunctionBegin;
121817ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1219416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
12203a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
1221e20fef11SSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
122217ab2063SBarry Smith #else
122317ab2063SBarry Smith       sum += (*v)*(*v); v++;
122417ab2063SBarry Smith #endif
122517ab2063SBarry Smith     }
122617ab2063SBarry Smith     *norm = sqrt(sum);
12273a40ed3dSBarry Smith   } else if (type == NORM_1) {
122817ab2063SBarry Smith     double *tmp;
1229416022c9SBarry Smith     int    *jj = a->j;
123066963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1231cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
123217ab2063SBarry Smith     *norm = 0.0;
1233416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1234a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
123517ab2063SBarry Smith     }
1236416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
123717ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
123817ab2063SBarry Smith     }
12390452661fSBarry Smith     PetscFree(tmp);
12403a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
124117ab2063SBarry Smith     *norm = 0.0;
1242416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1243416022c9SBarry Smith       v = a->a + a->i[j] + shift;
124417ab2063SBarry Smith       sum = 0.0;
1245416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1246cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
124717ab2063SBarry Smith       }
124817ab2063SBarry Smith       if (sum > *norm) *norm = sum;
124917ab2063SBarry Smith     }
12503a40ed3dSBarry Smith   } else {
1251a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
125217ab2063SBarry Smith   }
12533a40ed3dSBarry Smith   PetscFunctionReturn(0);
125417ab2063SBarry Smith }
125517ab2063SBarry Smith 
12565615d1e5SSatish Balay #undef __FUNC__
12575615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12588f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
125917ab2063SBarry Smith {
1260416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1261416022c9SBarry Smith   Mat        C;
1262416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1263416022c9SBarry Smith   int        shift = a->indexshift;
1264d5d45c9bSBarry Smith   Scalar     *array = a->a;
126517ab2063SBarry Smith 
12663a40ed3dSBarry Smith   PetscFunctionBegin;
1267a8c6a408SBarry Smith   if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place");
12680452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1269cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
127017ab2063SBarry Smith   if (shift) {
127117ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
127217ab2063SBarry Smith   }
127317ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1274416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12750452661fSBarry Smith   PetscFree(col);
127617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
127717ab2063SBarry Smith     len = ai[i+1]-ai[i];
1278416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
127917ab2063SBarry Smith     array += len; aj += len;
128017ab2063SBarry Smith   }
128117ab2063SBarry Smith   if (shift) {
128217ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
128317ab2063SBarry Smith   }
128417ab2063SBarry Smith 
12856d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
12866d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
128717ab2063SBarry Smith 
12883638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1289416022c9SBarry Smith     *B = C;
129017ab2063SBarry Smith   } else {
1291f830108cSBarry Smith     PetscOps *Abops;
12920a6ffc59SBarry Smith     MatOps   Aops;
1293f830108cSBarry Smith 
1294416022c9SBarry Smith     /* This isn't really an in-place transpose */
12950452661fSBarry Smith     PetscFree(a->a);
12960452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
12970452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
12980452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
12990452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
13000452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
13011073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
13020452661fSBarry Smith     PetscFree(a);
1303f830108cSBarry Smith 
1304488ecbafSBarry Smith 
1305488ecbafSBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
1306488ecbafSBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
1307488ecbafSBarry Smith 
1308f830108cSBarry Smith     /*
1309f830108cSBarry Smith       This is horrible, horrible code. We need to keep the
13108f0f457cSSatish Balay       the bops and ops Structures, copy everything from C
13118f0f457cSSatish Balay       including the function pointers..
1312f830108cSBarry Smith     */
13138f0f457cSSatish Balay     PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps));
13148f0f457cSSatish Balay     PetscMemcpy(A->bops,C->bops,sizeof(PetscOps));
1315f830108cSBarry Smith     Abops    = A->bops;
1316f830108cSBarry Smith     Aops     = A->ops;
1317f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
1318f830108cSBarry Smith     A->bops  = Abops;
1319f830108cSBarry Smith     A->ops   = Aops;
132027a8da17SBarry Smith     A->qlist = 0;
1321f830108cSBarry Smith 
13220452661fSBarry Smith     PetscHeaderDestroy(C);
132317ab2063SBarry Smith   }
13243a40ed3dSBarry Smith   PetscFunctionReturn(0);
132517ab2063SBarry Smith }
132617ab2063SBarry Smith 
13275615d1e5SSatish Balay #undef __FUNC__
13285615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
13298f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
133017ab2063SBarry Smith {
1331416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
133217ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1333e1311b90SBarry Smith   int        ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
133417ab2063SBarry Smith 
13353a40ed3dSBarry Smith   PetscFunctionBegin;
133617ab2063SBarry Smith   if (ll) {
13373ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
13383ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
1339e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr);
1340a8c6a408SBarry Smith     if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
1341e1311b90SBarry Smith     ierr = VecGetArray(ll,&l); CHKERRQ(ierr);
1342416022c9SBarry Smith     v = a->a;
134317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
134417ab2063SBarry Smith       x = l[i];
1345416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
134617ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
134717ab2063SBarry Smith     }
1348e1311b90SBarry Smith     ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr);
134944cd7ae7SLois Curfman McInnes     PLogFlops(nz);
135017ab2063SBarry Smith   }
135117ab2063SBarry Smith   if (rr) {
1352e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr);
1353a8c6a408SBarry Smith     if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
1354e1311b90SBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1355416022c9SBarry Smith     v = a->a; jj = a->j;
135617ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
135717ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
135817ab2063SBarry Smith     }
1359e1311b90SBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
136044cd7ae7SLois Curfman McInnes     PLogFlops(nz);
136117ab2063SBarry Smith   }
13623a40ed3dSBarry Smith   PetscFunctionReturn(0);
136317ab2063SBarry Smith }
136417ab2063SBarry Smith 
13655615d1e5SSatish Balay #undef __FUNC__
13665615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
13676a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B)
136817ab2063SBarry Smith {
1369db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1370d5db1b72SBarry Smith   int          *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
137199141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1372a2744918SBarry Smith   register int sum,lensi;
137399141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
137499141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
137599141d43SSatish Balay   Scalar       *a_new,*mat_a;
1376416022c9SBarry Smith   Mat          C;
1377fee21e36SBarry Smith   PetscTruth   stride;
137817ab2063SBarry Smith 
13793a40ed3dSBarry Smith   PetscFunctionBegin;
1380d64ed03dSBarry Smith   ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr);
1381a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
1382d64ed03dSBarry Smith   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
1383a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
138499141d43SSatish Balay 
138517ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
138617ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
138717ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
138817ab2063SBarry Smith 
1389fee21e36SBarry Smith   ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr);
1390fee21e36SBarry Smith   ierr = ISStride(iscol,&stride); CHKERRQ(ierr);
1391fee21e36SBarry Smith   if (stride && step == 1) {
139202834360SBarry Smith     /* special case of contiguous rows */
139357faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
139402834360SBarry Smith     starts = lens + ncols;
139502834360SBarry Smith     /* loop over new rows determining lens and starting points */
139602834360SBarry Smith     for (i=0; i<nrows; i++) {
1397a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1398a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
139902834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1400d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
140102834360SBarry Smith           starts[i] = k;
140202834360SBarry Smith           break;
140302834360SBarry Smith 	}
140402834360SBarry Smith       }
1405a2744918SBarry Smith       sum = 0;
140602834360SBarry Smith       while (k < kend) {
1407d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1408a2744918SBarry Smith         sum++;
140902834360SBarry Smith       }
1410a2744918SBarry Smith       lens[i] = sum;
141102834360SBarry Smith     }
141202834360SBarry Smith     /* create submatrix */
1413cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
141408480c60SBarry Smith       int n_cols,n_rows;
141508480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1416a8c6a408SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1417d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
141808480c60SBarry Smith       C = *B;
14193a40ed3dSBarry Smith     } else {
142002834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
142108480c60SBarry Smith     }
1422db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1423db02288aSLois Curfman McInnes 
142402834360SBarry Smith     /* loop over rows inserting into submatrix */
1425db02288aSLois Curfman McInnes     a_new    = c->a;
1426db02288aSLois Curfman McInnes     j_new    = c->j;
1427db02288aSLois Curfman McInnes     i_new    = c->i;
1428db02288aSLois Curfman McInnes     i_new[0] = -shift;
142902834360SBarry Smith     for (i=0; i<nrows; i++) {
1430a2744918SBarry Smith       ii    = starts[i];
1431a2744918SBarry Smith       lensi = lens[i];
1432a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1433a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
143402834360SBarry Smith       }
1435a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1436a2744918SBarry Smith       a_new      += lensi;
1437a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1438a2744918SBarry Smith       c->ilen[i]  = lensi;
143902834360SBarry Smith     }
14400452661fSBarry Smith     PetscFree(lens);
14413a40ed3dSBarry Smith   } else {
144202834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
14430452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
144402834360SBarry Smith     ssmap = smap + shift;
144599141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1446cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
144717ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
144802834360SBarry Smith     /* determine lens of each row */
144902834360SBarry Smith     for (i=0; i<nrows; i++) {
1450d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
145102834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
145202834360SBarry Smith       lens[i] = 0;
145302834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1454d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
145502834360SBarry Smith           lens[i]++;
145602834360SBarry Smith         }
145702834360SBarry Smith       }
145802834360SBarry Smith     }
145917ab2063SBarry Smith     /* Create and fill new matrix */
1460a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
146199141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
1462a8c6a408SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
146399141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1464a8c6a408SBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
146599141d43SSatish Balay       }
146699141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
146708480c60SBarry Smith       C = *B;
14683a40ed3dSBarry Smith     } else {
146902834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
147008480c60SBarry Smith     }
147199141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
147217ab2063SBarry Smith     for (i=0; i<nrows; i++) {
147399141d43SSatish Balay       row    = irow[i];
147499141d43SSatish Balay       kstart = ai[row]+shift;
147599141d43SSatish Balay       kend   = kstart + a->ilen[row];
147699141d43SSatish Balay       mat_i  = c->i[i]+shift;
147799141d43SSatish Balay       mat_j  = c->j + mat_i;
147899141d43SSatish Balay       mat_a  = c->a + mat_i;
147999141d43SSatish Balay       mat_ilen = c->ilen + i;
148017ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
148199141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
148299141d43SSatish Balay           *mat_j++ = tcol - (!shift);
148399141d43SSatish Balay           *mat_a++ = a->a[k];
148499141d43SSatish Balay           (*mat_ilen)++;
148599141d43SSatish Balay 
148617ab2063SBarry Smith         }
148717ab2063SBarry Smith       }
148817ab2063SBarry Smith     }
148902834360SBarry Smith     /* Free work space */
149002834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
149199141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
149202834360SBarry Smith   }
14936d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
14946d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
149517ab2063SBarry Smith 
149617ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1497416022c9SBarry Smith   *B = C;
14983a40ed3dSBarry Smith   PetscFunctionReturn(0);
149917ab2063SBarry Smith }
150017ab2063SBarry Smith 
1501a871dcd8SBarry Smith /*
150263b91edcSBarry Smith      note: This can only work for identity for row and col. It would
150363b91edcSBarry Smith    be good to check this and otherwise generate an error.
1504a871dcd8SBarry Smith */
15055615d1e5SSatish Balay #undef __FUNC__
15065615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
15078f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1508a871dcd8SBarry Smith {
150963b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
151008480c60SBarry Smith   int        ierr;
151163b91edcSBarry Smith   Mat        outA;
151263b91edcSBarry Smith 
15133a40ed3dSBarry Smith   PetscFunctionBegin;
1514a8c6a408SBarry Smith   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
1515a871dcd8SBarry Smith 
151663b91edcSBarry Smith   outA          = inA;
151763b91edcSBarry Smith   inA->factor   = FACTOR_LU;
151863b91edcSBarry Smith   a->row        = row;
151963b91edcSBarry Smith   a->col        = col;
152063b91edcSBarry Smith 
1521f0ec6fceSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1522f0ec6fceSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
15231e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1524f0ec6fceSSatish Balay 
152594a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
15260452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work);
152794a9d846SBarry Smith   }
152863b91edcSBarry Smith 
152908480c60SBarry Smith   if (!a->diag) {
153008480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
153163b91edcSBarry Smith   }
153263b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
15333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1534a871dcd8SBarry Smith }
1535a871dcd8SBarry Smith 
153674cdf0dfSBarry Smith #include "pinclude/blaslapack.h"
15375615d1e5SSatish Balay #undef __FUNC__
15385615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
15398f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1540f0b747eeSBarry Smith {
1541f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1542f0b747eeSBarry Smith   int        one = 1;
15433a40ed3dSBarry Smith 
15443a40ed3dSBarry Smith   PetscFunctionBegin;
1545f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1546f0b747eeSBarry Smith   PLogFlops(a->nz);
15473a40ed3dSBarry Smith   PetscFunctionReturn(0);
1548f0b747eeSBarry Smith }
1549f0b747eeSBarry Smith 
15505615d1e5SSatish Balay #undef __FUNC__
15515615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
15526a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1553cddf8d76SBarry Smith {
1554cddf8d76SBarry Smith   int ierr,i;
1555cddf8d76SBarry Smith 
15563a40ed3dSBarry Smith   PetscFunctionBegin;
1557cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
15580452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1559cddf8d76SBarry Smith   }
1560cddf8d76SBarry Smith 
1561cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
15626a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1563cddf8d76SBarry Smith   }
15643a40ed3dSBarry Smith   PetscFunctionReturn(0);
1565cddf8d76SBarry Smith }
1566cddf8d76SBarry Smith 
15675615d1e5SSatish Balay #undef __FUNC__
15685615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
15698f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
15705a838052SSatish Balay {
1571f830108cSBarry Smith   PetscFunctionBegin;
15725a838052SSatish Balay   *bs = 1;
15733a40ed3dSBarry Smith   PetscFunctionReturn(0);
15745a838052SSatish Balay }
15755a838052SSatish Balay 
15765615d1e5SSatish Balay #undef __FUNC__
15775615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15788f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15794dcbc457SBarry Smith {
1580e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
158106763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15828a047759SSatish Balay   int        start, end, *ai, *aj;
1583bbd702dbSSatish Balay   BT         table;
1584bbd702dbSSatish Balay 
15853a40ed3dSBarry Smith   PetscFunctionBegin;
15868a047759SSatish Balay   shift = a->indexshift;
1587e4d965acSSatish Balay   m     = a->m;
1588e4d965acSSatish Balay   ai    = a->i;
15898a047759SSatish Balay   aj    = a->j+shift;
15908a047759SSatish Balay 
1591a8c6a408SBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used");
159206763907SSatish Balay 
159306763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1594bbd702dbSSatish Balay   ierr  = BTCreate(m,table); CHKERRQ(ierr);
159506763907SSatish Balay 
1596e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1597b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1598e4d965acSSatish Balay     isz  = 0;
1599bbd702dbSSatish Balay     BTMemzero(m,table);
1600e4d965acSSatish Balay 
1601e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
16024dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
160377c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1604e4d965acSSatish Balay 
1605dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1606e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
1607bbd702dbSSatish Balay       if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];}
16084dcbc457SBarry Smith     }
160906763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
161006763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1611e4d965acSSatish Balay 
161204a348a9SBarry Smith     k = 0;
161304a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
161404a348a9SBarry Smith       n = isz;
161506763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1616e4d965acSSatish Balay         row   = nidx[k];
1617e4d965acSSatish Balay         start = ai[row];
1618e4d965acSSatish Balay         end   = ai[row+1];
161904a348a9SBarry Smith         for ( l = start; l<end ; l++){
16208a047759SSatish Balay           val = aj[l] + shift;
16212a8f2162SSatish Balay           if (!BTLookupSet(table,val)) {nidx[isz++] = val;}
1622e4d965acSSatish Balay         }
1623e4d965acSSatish Balay       }
1624e4d965acSSatish Balay     }
1625029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1626e4d965acSSatish Balay   }
1627bbd702dbSSatish Balay   BTDestroy(table);
162806763907SSatish Balay   PetscFree(nidx);
16293a40ed3dSBarry Smith   PetscFunctionReturn(0);
16304dcbc457SBarry Smith }
163117ab2063SBarry Smith 
16320513a670SBarry Smith /* -------------------------------------------------------------- */
16330513a670SBarry Smith #undef __FUNC__
16340513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
16350513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
16360513a670SBarry Smith {
16370513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
16380513a670SBarry Smith   Scalar     *vwork;
16390513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
16400513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
164156cd22aeSBarry Smith   IS         icolp,irowp;
16420513a670SBarry Smith 
16433a40ed3dSBarry Smith   PetscFunctionBegin;
164456cd22aeSBarry Smith   ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr);
164556cd22aeSBarry Smith   ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr);
164656cd22aeSBarry Smith   ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr);
164756cd22aeSBarry Smith   ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr);
16480513a670SBarry Smith 
16490513a670SBarry Smith   /* determine lengths of permuted rows */
16500513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
16510513a670SBarry Smith   for (i=0; i<m; i++ ) {
16520513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
16530513a670SBarry Smith   }
16540513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
16550513a670SBarry Smith   PetscFree(lens);
16560513a670SBarry Smith 
16570513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
16580513a670SBarry Smith   for (i=0; i<m; i++) {
16590513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16600513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
16610513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
16620513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16630513a670SBarry Smith   }
16640513a670SBarry Smith   PetscFree(cnew);
16650513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
16660513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
166756cd22aeSBarry Smith   ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr);
166856cd22aeSBarry Smith   ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr);
166956cd22aeSBarry Smith   ierr = ISDestroy(irowp); CHKERRQ(ierr);
167056cd22aeSBarry Smith   ierr = ISDestroy(icolp); CHKERRQ(ierr);
16713a40ed3dSBarry Smith   PetscFunctionReturn(0);
16720513a670SBarry Smith }
16730513a670SBarry Smith 
16745615d1e5SSatish Balay #undef __FUNC__
16755615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1676682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1677682d7d0cSBarry Smith {
1678682d7d0cSBarry Smith   static int called = 0;
1679682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1680682d7d0cSBarry Smith 
16813a40ed3dSBarry Smith   PetscFunctionBegin;
16823a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
168376be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
168476be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
168576be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
168676be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodes\n");
168776be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1688682d7d0cSBarry Smith #if defined(HAVE_ESSL)
168976be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1690682d7d0cSBarry Smith #endif
16913a40ed3dSBarry Smith   PetscFunctionReturn(0);
1692682d7d0cSBarry Smith }
16938f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1694a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1695a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1696a93ec695SBarry Smith 
1697*cb5b572fSBarry Smith #undef __FUNC__
1698*cb5b572fSBarry Smith #define __FUNC__ "MatCopy_SeqAIJ"
1699*cb5b572fSBarry Smith int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1700*cb5b572fSBarry Smith {
1701*cb5b572fSBarry Smith   int    ierr,i,rstart,rend,nz,*cwork;
1702*cb5b572fSBarry Smith   Scalar *vwork;
1703*cb5b572fSBarry Smith 
1704*cb5b572fSBarry Smith   PetscFunctionBegin;
1705*cb5b572fSBarry Smith   if (str == SAME_NONZERO_PATTERN) {
1706*cb5b572fSBarry Smith     ierr = MatZeroEntries(B); CHKERRQ(ierr);
1707*cb5b572fSBarry Smith     ierr = MatGetOwnershipRange(A,&rstart,&rend); CHKERRQ(ierr);
1708*cb5b572fSBarry Smith     for (i=rstart; i<rend; i++) {
1709*cb5b572fSBarry Smith       ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
1710*cb5b572fSBarry Smith       ierr = MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
1711*cb5b572fSBarry Smith       ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
1712*cb5b572fSBarry Smith     }
1713*cb5b572fSBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1714*cb5b572fSBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1715*cb5b572fSBarry Smith   } else {
1716*cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1717*cb5b572fSBarry Smith   }
1718*cb5b572fSBarry Smith   PetscFunctionReturn(0);
1719*cb5b572fSBarry Smith }
1720*cb5b572fSBarry Smith 
1721682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
17220a6ffc59SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
1723*cb5b572fSBarry Smith        MatGetRow_SeqAIJ,
1724*cb5b572fSBarry Smith        MatRestoreRow_SeqAIJ,
1725*cb5b572fSBarry Smith        MatMult_SeqAIJ,
1726*cb5b572fSBarry Smith        MatMultAdd_SeqAIJ,
1727*cb5b572fSBarry Smith        MatMultTrans_SeqAIJ,
1728*cb5b572fSBarry Smith        MatMultTransAdd_SeqAIJ,
1729*cb5b572fSBarry Smith        MatSolve_SeqAIJ,
1730*cb5b572fSBarry Smith        MatSolveAdd_SeqAIJ,
1731*cb5b572fSBarry Smith        MatSolveTrans_SeqAIJ,
1732*cb5b572fSBarry Smith        MatSolveTransAdd_SeqAIJ,
1733*cb5b572fSBarry Smith        MatLUFactor_SeqAIJ,
1734*cb5b572fSBarry Smith        0,
173517ab2063SBarry Smith        MatRelax_SeqAIJ,
173617ab2063SBarry Smith        MatTranspose_SeqAIJ,
1737*cb5b572fSBarry Smith        MatGetInfo_SeqAIJ,
1738*cb5b572fSBarry Smith        MatEqual_SeqAIJ,
1739*cb5b572fSBarry Smith        MatGetDiagonal_SeqAIJ,
1740*cb5b572fSBarry Smith        MatDiagonalScale_SeqAIJ,
1741*cb5b572fSBarry Smith        MatNorm_SeqAIJ,
1742*cb5b572fSBarry Smith        0,
1743*cb5b572fSBarry Smith        MatAssemblyEnd_SeqAIJ,
174417ab2063SBarry Smith        MatCompress_SeqAIJ,
1745*cb5b572fSBarry Smith        MatSetOption_SeqAIJ,
1746*cb5b572fSBarry Smith        MatZeroEntries_SeqAIJ,
1747*cb5b572fSBarry Smith        MatZeroRows_SeqAIJ,
1748*cb5b572fSBarry Smith        MatLUFactorSymbolic_SeqAIJ,
1749*cb5b572fSBarry Smith        MatLUFactorNumeric_SeqAIJ,
1750*cb5b572fSBarry Smith        0,
1751*cb5b572fSBarry Smith        0,
1752*cb5b572fSBarry Smith        MatGetSize_SeqAIJ,
1753*cb5b572fSBarry Smith        MatGetSize_SeqAIJ,
1754*cb5b572fSBarry Smith        MatGetOwnershipRange_SeqAIJ,
1755*cb5b572fSBarry Smith        MatILUFactorSymbolic_SeqAIJ,
1756*cb5b572fSBarry Smith        0,
1757*cb5b572fSBarry Smith        0,
1758*cb5b572fSBarry Smith        0,
1759*cb5b572fSBarry Smith        MatConvertSameType_SeqAIJ,
1760*cb5b572fSBarry Smith        0,
1761*cb5b572fSBarry Smith        0,
1762*cb5b572fSBarry Smith        MatILUFactor_SeqAIJ,
1763*cb5b572fSBarry Smith        0,
1764*cb5b572fSBarry Smith        0,
1765*cb5b572fSBarry Smith        MatGetSubMatrices_SeqAIJ,
1766*cb5b572fSBarry Smith        MatIncreaseOverlap_SeqAIJ,
1767*cb5b572fSBarry Smith        MatGetValues_SeqAIJ,
1768*cb5b572fSBarry Smith        MatCopy_SeqAIJ,
1769f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
1770*cb5b572fSBarry Smith        MatScale_SeqAIJ,
1771*cb5b572fSBarry Smith        0,
1772*cb5b572fSBarry Smith        0,
17736945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
17746945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
17753b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
17763b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
17773b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1778a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1779a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
17800513a670SBarry Smith        MatColoringPatch_SeqAIJ,
17810513a670SBarry Smith        0,
1782cda55fadSBarry Smith        MatPermute_SeqAIJ,
1783cda55fadSBarry Smith        0,
1784cda55fadSBarry Smith        0,
1785cda55fadSBarry Smith        0,
1786cda55fadSBarry Smith        0,
1787cda55fadSBarry Smith        MatGetMaps_Petsc};
178817ab2063SBarry Smith 
178917ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
179017ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
179117ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
179217ab2063SBarry Smith 
17935615d1e5SSatish Balay #undef __FUNC__
1794bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ"
1795bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
1796bef8e0ddSBarry Smith {
1797bef8e0ddSBarry Smith   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
1798bef8e0ddSBarry Smith   int        i,nz,n;
1799bef8e0ddSBarry Smith 
1800bef8e0ddSBarry Smith   PetscFunctionBegin;
1801bef8e0ddSBarry Smith   if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing");
1802bef8e0ddSBarry Smith 
1803bef8e0ddSBarry Smith   nz = aij->maxnz;
1804bef8e0ddSBarry Smith   n  = aij->n;
1805bef8e0ddSBarry Smith   for (i=0; i<nz; i++) {
1806bef8e0ddSBarry Smith     aij->j[i] = indices[i];
1807bef8e0ddSBarry Smith   }
1808bef8e0ddSBarry Smith   aij->nz = nz;
1809bef8e0ddSBarry Smith   for ( i=0; i<n; i++ ) {
1810bef8e0ddSBarry Smith     aij->ilen[i] = aij->imax[i];
1811bef8e0ddSBarry Smith   }
1812bef8e0ddSBarry Smith 
1813bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1814bef8e0ddSBarry Smith }
1815bef8e0ddSBarry Smith 
1816bef8e0ddSBarry Smith #undef __FUNC__
1817bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices"
1818bef8e0ddSBarry Smith /*@
1819bef8e0ddSBarry Smith     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
1820bef8e0ddSBarry Smith        in the matrix.
1821bef8e0ddSBarry Smith 
1822bef8e0ddSBarry Smith   Input Parameters:
1823bef8e0ddSBarry Smith +  mat - the SeqAIJ matrix
1824bef8e0ddSBarry Smith -  indices - the column indices
1825bef8e0ddSBarry Smith 
1826bef8e0ddSBarry Smith   Notes:
1827bef8e0ddSBarry Smith     This can be called if you have precomputed the nonzero structure of the
1828bef8e0ddSBarry Smith   matrix and want to provide it to the matrix object to improve the performance
1829bef8e0ddSBarry Smith   of the MatSetValues() operation.
1830bef8e0ddSBarry Smith 
1831bef8e0ddSBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
1832bef8e0ddSBarry Smith   MatCreateSeqAIJ().
1833bef8e0ddSBarry Smith 
1834bef8e0ddSBarry Smith     MUST be called before any calls to MatSetValues();
1835bef8e0ddSBarry Smith 
1836bef8e0ddSBarry Smith @*/
1837bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
1838bef8e0ddSBarry Smith {
1839bef8e0ddSBarry Smith   int ierr,(*f)(Mat,int *);
1840bef8e0ddSBarry Smith 
1841bef8e0ddSBarry Smith   PetscFunctionBegin;
1842bef8e0ddSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1843bef8e0ddSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f);
1844bef8e0ddSBarry Smith          CHKERRQ(ierr);
1845bef8e0ddSBarry Smith   if (f) {
1846bef8e0ddSBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
1847bef8e0ddSBarry Smith   } else {
1848bef8e0ddSBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
1849bef8e0ddSBarry Smith   }
1850bef8e0ddSBarry Smith   PetscFunctionReturn(0);
1851bef8e0ddSBarry Smith }
1852bef8e0ddSBarry Smith 
1853*cb5b572fSBarry Smith 
1854bef8e0ddSBarry Smith #undef __FUNC__
18555615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
185617ab2063SBarry Smith /*@C
1857682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
18580d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
18596e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
18602bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
18612bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
186217ab2063SBarry Smith 
1863db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1864db81eaa0SLois Curfman McInnes 
186517ab2063SBarry Smith    Input Parameters:
1866db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
186717ab2063SBarry Smith .  m - number of rows
186817ab2063SBarry Smith .  n - number of columns
186917ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
1870db81eaa0SLois Curfman McInnes -  nzz - array containing the number of nonzeros in the various rows
18712bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
187217ab2063SBarry Smith 
187317ab2063SBarry Smith    Output Parameter:
1874416022c9SBarry Smith .  A - the matrix
187517ab2063SBarry Smith 
1876b259b22eSLois Curfman McInnes    Notes:
187717ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
187817ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
18790002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
188044cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
188117ab2063SBarry Smith 
188217ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1883a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
18843d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
18856da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
188617ab2063SBarry Smith 
1887682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
18884fca80b9SLois Curfman McInnes    improve numerical efficiency of matrix-vector products and solves. We
1889682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
18906c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
18916c7ebb05SLois Curfman McInnes 
18926c7ebb05SLois Curfman McInnes    Options Database Keys:
1893db81eaa0SLois Curfman McInnes +  -mat_aij_no_inode  - Do not use inodes
1894db81eaa0SLois Curfman McInnes .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
1895db81eaa0SLois Curfman McInnes -  -mat_aij_oneindex - Internally use indexing starting at 1
1896db81eaa0SLois Curfman McInnes         rather than 0.  Note that when calling MatSetValues(),
1897db81eaa0SLois Curfman McInnes         the user still MUST index entries starting at 0!
189817ab2063SBarry Smith 
1899bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices()
190017ab2063SBarry Smith @*/
1901416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
190217ab2063SBarry Smith {
1903416022c9SBarry Smith   Mat        B;
1904416022c9SBarry Smith   Mat_SeqAIJ *b;
19056945ee14SBarry Smith   int        i, len, ierr, flg,size;
19066945ee14SBarry Smith 
19073a40ed3dSBarry Smith   PetscFunctionBegin;
19086945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1909a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1");
1910d5d45c9bSBarry Smith 
1911416022c9SBarry Smith   *A                  = 0;
1912f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView);
1913416022c9SBarry Smith   PLogObjectCreate(B);
19140452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
191544cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
19160a6ffc59SBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1917e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_SeqAIJ;
1918e1311b90SBarry Smith   B->ops->view             = MatView_SeqAIJ;
1919416022c9SBarry Smith   B->factor           = 0;
1920416022c9SBarry Smith   B->lupivotthreshold = 1.0;
192190f02eecSBarry Smith   B->mapping          = 0;
1922e1311b90SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg);CHKERRQ(ierr);
19237a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
1924e1311b90SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*)&b->ilu_preserve_row_sums);CHKERRQ(ierr);
1925416022c9SBarry Smith   b->row              = 0;
1926416022c9SBarry Smith   b->col              = 0;
192782bf6240SBarry Smith   b->icol             = 0;
1928416022c9SBarry Smith   b->indexshift       = 0;
1929b810aeb4SBarry Smith   b->reallocs         = 0;
193069957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
193169957df2SSatish Balay   if (flg) b->indexshift = -1;
193217ab2063SBarry Smith 
193344cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
193444cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
1935a5ae1ecdSBarry Smith 
19364d197716SBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
19374d197716SBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1938a5ae1ecdSBarry Smith 
19390452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1940b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
19417b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
19427b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1943416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
194417ab2063SBarry Smith     nz = nz*m;
19453a40ed3dSBarry Smith   } else {
194617ab2063SBarry Smith     nz = 0;
1947416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
194817ab2063SBarry Smith   }
194917ab2063SBarry Smith 
195017ab2063SBarry Smith   /* allocate the matrix space */
1951416022c9SBarry Smith   len             = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
19520452661fSBarry Smith   b->a            = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1953416022c9SBarry Smith   b->j            = (int *) (b->a + nz);
1954cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1955416022c9SBarry Smith   b->i            = b->j + nz;
1956416022c9SBarry Smith   b->singlemalloc = 1;
195717ab2063SBarry Smith 
1958416022c9SBarry Smith   b->i[0] = -b->indexshift;
195917ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1960416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
196117ab2063SBarry Smith   }
196217ab2063SBarry Smith 
1963416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
19640452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1965f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1966416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
196717ab2063SBarry Smith 
1968416022c9SBarry Smith   b->nz               = 0;
1969416022c9SBarry Smith   b->maxnz            = nz;
1970416022c9SBarry Smith   b->sorted           = 0;
1971416022c9SBarry Smith   b->roworiented      = 1;
1972416022c9SBarry Smith   b->nonew            = 0;
1973416022c9SBarry Smith   b->diag             = 0;
1974416022c9SBarry Smith   b->solve_work       = 0;
197571bd300dSLois Curfman McInnes   b->spptr            = 0;
1976754ec7b1SSatish Balay   b->inode.node_count = 0;
1977754ec7b1SSatish Balay   b->inode.size       = 0;
19786c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
19796c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
19804e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
198117ab2063SBarry Smith 
1982416022c9SBarry Smith   *A = B;
19834e220ebcSLois Curfman McInnes 
19844b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
19854b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
198669957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
198769957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
19884b14c69eSBarry Smith #endif
198969957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
199069957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
199169957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
199269957df2SSatish Balay   if (flg) {
1993a8c6a408SBarry Smith     if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1994416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
199517ab2063SBarry Smith   }
199669957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
199769957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1998bef8e0ddSBarry Smith 
1999bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2000bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2001bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
20023a40ed3dSBarry Smith   PetscFunctionReturn(0);
200317ab2063SBarry Smith }
200417ab2063SBarry Smith 
20055615d1e5SSatish Balay #undef __FUNC__
20065615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
20073d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
200817ab2063SBarry Smith {
2009416022c9SBarry Smith   Mat        C;
2010416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
2011bef8e0ddSBarry Smith   int        i,len, m = a->m,shift = a->indexshift,ierr;
201217ab2063SBarry Smith 
20133a40ed3dSBarry Smith   PetscFunctionBegin;
20144043dd9cSLois Curfman McInnes   *B = 0;
2015f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView);
2016416022c9SBarry Smith   PLogObjectCreate(C);
20170452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
2018f830108cSBarry Smith   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2019e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqAIJ;
2020e1311b90SBarry Smith   C->ops->view       = MatView_SeqAIJ;
2021416022c9SBarry Smith   C->factor     = A->factor;
2022416022c9SBarry Smith   c->row        = 0;
2023416022c9SBarry Smith   c->col        = 0;
202482bf6240SBarry Smith   c->icol       = 0;
2025416022c9SBarry Smith   c->indexshift = shift;
2026c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
202717ab2063SBarry Smith 
202844cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
202944cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
203044cd7ae7SLois Curfman McInnes   C->M          = a->m;
203144cd7ae7SLois Curfman McInnes   C->N          = a->n;
203217ab2063SBarry Smith 
20330452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
20340452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
203517ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
2036416022c9SBarry Smith     c->imax[i] = a->imax[i];
2037416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
203817ab2063SBarry Smith   }
203917ab2063SBarry Smith 
204017ab2063SBarry Smith   /* allocate the matrix space */
2041416022c9SBarry Smith   c->singlemalloc = 1;
2042416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
20430452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
2044416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
2045416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
2046416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
204717ab2063SBarry Smith   if (m > 0) {
2048416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
204908480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
2050416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
205117ab2063SBarry Smith     }
205208480c60SBarry Smith   }
205317ab2063SBarry Smith 
2054f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2055416022c9SBarry Smith   c->sorted      = a->sorted;
2056416022c9SBarry Smith   c->roworiented = a->roworiented;
2057416022c9SBarry Smith   c->nonew       = a->nonew;
20587a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
205917ab2063SBarry Smith 
2060416022c9SBarry Smith   if (a->diag) {
20610452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
2062416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
206317ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
2064416022c9SBarry Smith       c->diag[i] = a->diag[i];
206517ab2063SBarry Smith     }
20663a40ed3dSBarry Smith   } else c->diag          = 0;
20676c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
20686c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
2069754ec7b1SSatish Balay   if (a->inode.size){
2070daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
2071754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
2072daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
2073754ec7b1SSatish Balay   } else {
2074754ec7b1SSatish Balay     c->inode.size       = 0;
2075754ec7b1SSatish Balay     c->inode.node_count = 0;
2076754ec7b1SSatish Balay   }
2077416022c9SBarry Smith   c->nz                 = a->nz;
2078416022c9SBarry Smith   c->maxnz              = a->maxnz;
2079416022c9SBarry Smith   c->solve_work         = 0;
208076dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
2081754ec7b1SSatish Balay 
2082416022c9SBarry Smith   *B = C;
2083bef8e0ddSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqAIJSetColumnIndices_C",
2084bef8e0ddSBarry Smith                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2085bef8e0ddSBarry Smith                                      (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr);
20863a40ed3dSBarry Smith   PetscFunctionReturn(0);
208717ab2063SBarry Smith }
208817ab2063SBarry Smith 
20895615d1e5SSatish Balay #undef __FUNC__
20905615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
209119bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
209217ab2063SBarry Smith {
2093416022c9SBarry Smith   Mat_SeqAIJ   *a;
2094416022c9SBarry Smith   Mat          B;
209517699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
2096bcd2baecSBarry Smith   MPI_Comm     comm;
209717ab2063SBarry Smith 
20983a40ed3dSBarry Smith   PetscFunctionBegin;
209919bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
210017699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
2101a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor");
210290ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
21030752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
2104a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file");
210517ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
210617ab2063SBarry Smith 
2107d64ed03dSBarry Smith   if (nz < 0) {
2108a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ");
2109d64ed03dSBarry Smith   }
2110d64ed03dSBarry Smith 
211117ab2063SBarry Smith   /* read in row lengths */
21120452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
21130752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
211417ab2063SBarry Smith 
211517ab2063SBarry Smith   /* create our matrix */
2116416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
2117416022c9SBarry Smith   B = *A;
2118416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
2119416022c9SBarry Smith   shift = a->indexshift;
212017ab2063SBarry Smith 
212117ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
21220752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr);
212317ab2063SBarry Smith   if (shift) {
212417ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
2125416022c9SBarry Smith       a->j[i] += 1;
212617ab2063SBarry Smith     }
212717ab2063SBarry Smith   }
212817ab2063SBarry Smith 
212917ab2063SBarry Smith   /* read in nonzero values */
21300752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr);
213117ab2063SBarry Smith 
213217ab2063SBarry Smith   /* set matrix "i" values */
2133416022c9SBarry Smith   a->i[0] = -shift;
213417ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
2135416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2136416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
213717ab2063SBarry Smith   }
21380452661fSBarry Smith   PetscFree(rowlengths);
213917ab2063SBarry Smith 
21406d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
21416d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
21423a40ed3dSBarry Smith   PetscFunctionReturn(0);
214317ab2063SBarry Smith }
214417ab2063SBarry Smith 
21455615d1e5SSatish Balay #undef __FUNC__
21465615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
21478f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
21487264ac53SSatish Balay {
21497264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
21507264ac53SSatish Balay 
21513a40ed3dSBarry Smith   PetscFunctionBegin;
2152a8c6a408SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
21537264ac53SSatish Balay 
21547264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
21557264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
2156bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
21573a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2158bcd2baecSBarry Smith   }
21597264ac53SSatish Balay 
21607264ac53SSatish Balay   /* if the a->i are the same */
21618108c231SLois Curfman McInnes   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
21623a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
21637264ac53SSatish Balay   }
21647264ac53SSatish Balay 
21657264ac53SSatish Balay   /* if a->j are the same */
2166bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
21673a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2168bcd2baecSBarry Smith   }
2169bcd2baecSBarry Smith 
2170bcd2baecSBarry Smith   /* if a->a are the same */
217119bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
21723a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
21737264ac53SSatish Balay   }
217477c4ece6SBarry Smith   *flg = PETSC_TRUE;
21753a40ed3dSBarry Smith   PetscFunctionReturn(0);
21767264ac53SSatish Balay 
21777264ac53SSatish Balay }
2178