xref: /petsc/src/mat/impls/aij/seq/aij.c (revision e1311b9049e89cb3452dcd306fde571f4b440ff2)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*e1311b90SBarry Smith static char vcid[] = "$Id: aij.c,v 1.255 1998/03/26 22:32:03 balay Exp bsmith $";
317ab2063SBarry Smith #endif
417ab2063SBarry Smith 
5d5d45c9bSBarry Smith /*
63369ce9aSBarry Smith     Defines the basic matrix operations for the AIJ (compressed row)
7d5d45c9bSBarry Smith   matrix storage format.
8d5d45c9bSBarry Smith */
93369ce9aSBarry Smith 
103369ce9aSBarry Smith #include "pinclude/pviewer.h"
113369ce9aSBarry Smith #include "sys.h"
1270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
13f5eb4b81SSatish Balay #include "src/vec/vecimpl.h"
14f5eb4b81SSatish Balay #include "src/inline/spops.h"
158d195f9aSBarry Smith #include "src/inline/dot.h"
16f5eb4b81SSatish Balay #include "src/inline/bitarray.h"
1717ab2063SBarry Smith 
18a2ce50c7SBarry Smith /*
19a2ce50c7SBarry Smith     Basic AIJ format ILU based on drop tolerance
20a2ce50c7SBarry Smith */
215615d1e5SSatish Balay #undef __FUNC__
225615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ"
23a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact)
24a2ce50c7SBarry Smith {
25a2ce50c7SBarry Smith   /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */
26a2ce50c7SBarry Smith   int        ierr = 1;
27a2ce50c7SBarry Smith 
283a40ed3dSBarry Smith   PetscFunctionBegin;
29e3372554SBarry Smith   SETERRQ(ierr,0,"Not implemented");
30d64ed03dSBarry Smith #if !defined(USE_PETSC_DEBUG)
31d64ed03dSBarry Smith   PetscFunctionReturn(0);
32d64ed03dSBarry Smith #endif
33a2ce50c7SBarry Smith }
34a2ce50c7SBarry Smith 
35bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
3617ab2063SBarry Smith 
375615d1e5SSatish Balay #undef __FUNC__
385615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ"
398f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,
406945ee14SBarry Smith                            PetscTruth *done)
4117ab2063SBarry Smith {
42416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
436945ee14SBarry Smith   int        ierr,i,ishift;
4417ab2063SBarry Smith 
453a40ed3dSBarry Smith   PetscFunctionBegin;
4631625ec5SSatish Balay   *m     = A->m;
473a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
486945ee14SBarry Smith   ishift = a->indexshift;
496945ee14SBarry Smith   if (symmetric) {
5031625ec5SSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
516945ee14SBarry Smith   } else if (oshift == 0 && ishift == -1) {
5231625ec5SSatish Balay     int nz = a->i[a->m];
533b2fbd54SBarry Smith     /* malloc space and  subtract 1 from i and j indices */
5431625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
553b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
563b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1;
5731625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1;
586945ee14SBarry Smith   } else if (oshift == 1 && ishift == 0) {
5931625ec5SSatish Balay     int nz = a->i[a->m] + 1;
603b2fbd54SBarry Smith     /* malloc space and  add 1 to i and j indices */
6131625ec5SSatish Balay     *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia);
623b2fbd54SBarry Smith     *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja);
633b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1;
6431625ec5SSatish Balay     for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1;
656945ee14SBarry Smith   } else {
666945ee14SBarry Smith     *ia = a->i; *ja = a->j;
67a2ce50c7SBarry Smith   }
68a2ce50c7SBarry Smith 
693a40ed3dSBarry Smith   PetscFunctionReturn(0);
70a2744918SBarry Smith }
71a2744918SBarry Smith 
725615d1e5SSatish Balay #undef __FUNC__
735615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ"
748f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,
756945ee14SBarry Smith                                PetscTruth *done)
766945ee14SBarry Smith {
776945ee14SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
783b2fbd54SBarry Smith   int        ishift = a->indexshift;
796945ee14SBarry Smith 
803a40ed3dSBarry Smith   PetscFunctionBegin;
813a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
823b2fbd54SBarry Smith   if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
836945ee14SBarry Smith     PetscFree(*ia);
846945ee14SBarry Smith     PetscFree(*ja);
85bcd2baecSBarry Smith   }
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
8717ab2063SBarry Smith }
8817ab2063SBarry Smith 
895615d1e5SSatish Balay #undef __FUNC__
905615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ"
9143a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
923b2fbd54SBarry Smith                            PetscTruth *done)
933b2fbd54SBarry Smith {
943b2fbd54SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
95a93ec695SBarry Smith   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
96a93ec695SBarry Smith   int        nz = a->i[m]+ishift,row,*jj,mr,col;
973b2fbd54SBarry Smith 
983a40ed3dSBarry Smith   PetscFunctionBegin;
993b2fbd54SBarry Smith   *nn     = A->n;
1003a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1013b2fbd54SBarry Smith   if (symmetric) {
102179192dfSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr);
1033b2fbd54SBarry Smith   } else {
10461d2ded1SBarry Smith     collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths);
1053b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1063b2fbd54SBarry Smith     cia        = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia);
107a93ec695SBarry Smith     cja        = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja);
1083b2fbd54SBarry Smith     jj = a->j;
1093b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) {
1103b2fbd54SBarry Smith       collengths[jj[i] + ishift]++;
1113b2fbd54SBarry Smith     }
1123b2fbd54SBarry Smith     cia[0] = oshift;
1133b2fbd54SBarry Smith     for ( i=0; i<n; i++) {
1143b2fbd54SBarry Smith       cia[i+1] = cia[i] + collengths[i];
1153b2fbd54SBarry Smith     }
1163b2fbd54SBarry Smith     PetscMemzero(collengths,n*sizeof(int));
1173b2fbd54SBarry Smith     jj = a->j;
118a93ec695SBarry Smith     for ( row=0; row<m; row++ ) {
119a93ec695SBarry Smith       mr = a->i[row+1] - a->i[row];
120a93ec695SBarry Smith       for ( i=0; i<mr; i++ ) {
1213b2fbd54SBarry Smith         col = *jj++ + ishift;
1223b2fbd54SBarry Smith         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1233b2fbd54SBarry Smith       }
1243b2fbd54SBarry Smith     }
1253b2fbd54SBarry Smith     PetscFree(collengths);
1263b2fbd54SBarry Smith     *ia = cia; *ja = cja;
1273b2fbd54SBarry Smith   }
1283b2fbd54SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1303b2fbd54SBarry Smith }
1313b2fbd54SBarry Smith 
1325615d1e5SSatish Balay #undef __FUNC__
1335615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ"
13443a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,
1353b2fbd54SBarry Smith                                      int **ja,PetscTruth *done)
1363b2fbd54SBarry Smith {
1373a40ed3dSBarry Smith   PetscFunctionBegin;
1383a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1393b2fbd54SBarry Smith 
1403b2fbd54SBarry Smith   PetscFree(*ia);
1413b2fbd54SBarry Smith   PetscFree(*ja);
1423b2fbd54SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1443b2fbd54SBarry Smith }
1453b2fbd54SBarry Smith 
146227d817aSBarry Smith #define CHUNKSIZE   15
14717ab2063SBarry Smith 
1485615d1e5SSatish Balay #undef __FUNC__
1495615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ"
15061d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
15117ab2063SBarry Smith {
152416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
153416022c9SBarry Smith   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
1544b0e389bSBarry Smith   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented;
155d5d45c9bSBarry Smith   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
156416022c9SBarry Smith   Scalar     *ap,value, *aa = a->a;
15717ab2063SBarry Smith 
1583a40ed3dSBarry Smith   PetscFunctionBegin;
15917ab2063SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
160416022c9SBarry Smith     row  = im[k];
1613a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
162a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
163a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
1643b2fbd54SBarry Smith #endif
16517ab2063SBarry Smith     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
16617ab2063SBarry Smith     rmax = imax[row]; nrow = ailen[row];
167416022c9SBarry Smith     low = 0;
16817ab2063SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
1693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
170a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
171a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
1723b2fbd54SBarry Smith #endif
1734b0e389bSBarry Smith       col = in[l] - shift;
1744b0e389bSBarry Smith       if (roworiented) {
1754b0e389bSBarry Smith         value = *v++;
1764b0e389bSBarry Smith       }
1774b0e389bSBarry Smith       else {
1784b0e389bSBarry Smith         value = v[k + l*m];
1794b0e389bSBarry Smith       }
180416022c9SBarry Smith       if (!sorted) low = 0; high = nrow;
181416022c9SBarry Smith       while (high-low > 5) {
182416022c9SBarry Smith         t = (low+high)/2;
183416022c9SBarry Smith         if (rp[t] > col) high = t;
184416022c9SBarry Smith         else             low  = t;
18517ab2063SBarry Smith       }
186416022c9SBarry Smith       for ( i=low; i<high; i++ ) {
18717ab2063SBarry Smith         if (rp[i] > col) break;
18817ab2063SBarry Smith         if (rp[i] == col) {
189416022c9SBarry Smith           if (is == ADD_VALUES) ap[i] += value;
19017ab2063SBarry Smith           else                  ap[i] = value;
19117ab2063SBarry Smith           goto noinsert;
19217ab2063SBarry Smith         }
19317ab2063SBarry Smith       }
194c2653b3dSLois Curfman McInnes       if (nonew == 1) goto noinsert;
195a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
19617ab2063SBarry Smith       if (nrow >= rmax) {
19717ab2063SBarry Smith         /* there is no extra room in row, therefore enlarge */
198416022c9SBarry Smith         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
19917ab2063SBarry Smith         Scalar *new_a;
20017ab2063SBarry Smith 
201a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
20296854ed6SLois Curfman McInnes 
20317ab2063SBarry Smith         /* malloc new storage space */
204416022c9SBarry Smith         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
2050452661fSBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
20617ab2063SBarry Smith         new_j   = (int *) (new_a + new_nz);
20717ab2063SBarry Smith         new_i   = new_j + new_nz;
20817ab2063SBarry Smith 
20917ab2063SBarry Smith         /* copy over old data into new slots */
21017ab2063SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
211416022c9SBarry Smith         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
212416022c9SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
213416022c9SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
214416022c9SBarry Smith         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
21517ab2063SBarry Smith                                                            len*sizeof(int));
216416022c9SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
217416022c9SBarry Smith         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
21817ab2063SBarry Smith                                                            len*sizeof(Scalar));
21917ab2063SBarry Smith         /* free up old matrix storage */
2200452661fSBarry Smith         PetscFree(a->a);
2210452661fSBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
222416022c9SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
223416022c9SBarry Smith         a->singlemalloc = 1;
22417ab2063SBarry Smith 
22517ab2063SBarry Smith         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
226416022c9SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
227416022c9SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
228416022c9SBarry Smith         a->maxnz += CHUNKSIZE;
229b810aeb4SBarry Smith         a->reallocs++;
23017ab2063SBarry Smith       }
231416022c9SBarry Smith       N = nrow++ - 1; a->nz++;
232416022c9SBarry Smith       /* shift up all the later entries in this row */
233416022c9SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
23417ab2063SBarry Smith         rp[ii+1] = rp[ii];
23517ab2063SBarry Smith         ap[ii+1] = ap[ii];
23617ab2063SBarry Smith       }
23717ab2063SBarry Smith       rp[i] = col;
23817ab2063SBarry Smith       ap[i] = value;
23917ab2063SBarry Smith       noinsert:;
240416022c9SBarry Smith       low = i + 1;
24117ab2063SBarry Smith     }
24217ab2063SBarry Smith     ailen[row] = nrow;
24317ab2063SBarry Smith   }
2443a40ed3dSBarry Smith   PetscFunctionReturn(0);
24517ab2063SBarry Smith }
24617ab2063SBarry Smith 
2475615d1e5SSatish Balay #undef __FUNC__
2485615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ"
2498f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
2507eb43aa7SLois Curfman McInnes {
2517eb43aa7SLois Curfman McInnes   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
252b49de8d1SLois Curfman McInnes   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
2537eb43aa7SLois Curfman McInnes   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
2547eb43aa7SLois Curfman McInnes   Scalar     *ap, *aa = a->a, zero = 0.0;
2557eb43aa7SLois Curfman McInnes 
2563a40ed3dSBarry Smith   PetscFunctionBegin;
2577eb43aa7SLois Curfman McInnes   for ( k=0; k<m; k++ ) { /* loop over rows */
2587eb43aa7SLois Curfman McInnes     row  = im[k];
259a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
260a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
2617eb43aa7SLois Curfman McInnes     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
2627eb43aa7SLois Curfman McInnes     nrow = ailen[row];
2637eb43aa7SLois Curfman McInnes     for ( l=0; l<n; l++ ) { /* loop over columns */
264a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
265a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
2667eb43aa7SLois Curfman McInnes       col = in[l] - shift;
2677eb43aa7SLois Curfman McInnes       high = nrow; low = 0; /* assume unsorted */
2687eb43aa7SLois Curfman McInnes       while (high-low > 5) {
2697eb43aa7SLois Curfman McInnes         t = (low+high)/2;
2707eb43aa7SLois Curfman McInnes         if (rp[t] > col) high = t;
2717eb43aa7SLois Curfman McInnes         else             low  = t;
2727eb43aa7SLois Curfman McInnes       }
2737eb43aa7SLois Curfman McInnes       for ( i=low; i<high; i++ ) {
2747eb43aa7SLois Curfman McInnes         if (rp[i] > col) break;
2757eb43aa7SLois Curfman McInnes         if (rp[i] == col) {
276b49de8d1SLois Curfman McInnes           *v++ = ap[i];
2777eb43aa7SLois Curfman McInnes           goto finished;
2787eb43aa7SLois Curfman McInnes         }
2797eb43aa7SLois Curfman McInnes       }
280b49de8d1SLois Curfman McInnes       *v++ = zero;
2817eb43aa7SLois Curfman McInnes       finished:;
2827eb43aa7SLois Curfman McInnes     }
2837eb43aa7SLois Curfman McInnes   }
2843a40ed3dSBarry Smith   PetscFunctionReturn(0);
2857eb43aa7SLois Curfman McInnes }
2867eb43aa7SLois Curfman McInnes 
28717ab2063SBarry Smith 
2885615d1e5SSatish Balay #undef __FUNC__
2895615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary"
2908f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
29117ab2063SBarry Smith {
292416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
293416022c9SBarry Smith   int        i, fd, *col_lens, ierr;
29417ab2063SBarry Smith 
2953a40ed3dSBarry Smith   PetscFunctionBegin;
29690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
2970452661fSBarry Smith   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
298416022c9SBarry Smith   col_lens[0] = MAT_COOKIE;
299416022c9SBarry Smith   col_lens[1] = a->m;
300416022c9SBarry Smith   col_lens[2] = a->n;
301416022c9SBarry Smith   col_lens[3] = a->nz;
302416022c9SBarry Smith 
303416022c9SBarry Smith   /* store lengths of each row and write (including header) to file */
304416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
305416022c9SBarry Smith     col_lens[4+i] = a->i[i+1] - a->i[i];
30617ab2063SBarry Smith   }
3070752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3080452661fSBarry Smith   PetscFree(col_lens);
309416022c9SBarry Smith 
310416022c9SBarry Smith   /* store column indices (zero start index) */
311416022c9SBarry Smith   if (a->indexshift) {
312416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]--;
31317ab2063SBarry Smith   }
3140752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0); CHKERRQ(ierr);
315416022c9SBarry Smith   if (a->indexshift) {
316416022c9SBarry Smith     for ( i=0; i<a->nz; i++ ) a->j[i]++;
31717ab2063SBarry Smith   }
318416022c9SBarry Smith 
319416022c9SBarry Smith   /* store nonzero values */
3200752156aSBarry Smith   ierr = PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
3213a40ed3dSBarry Smith   PetscFunctionReturn(0);
32217ab2063SBarry Smith }
323416022c9SBarry Smith 
3245615d1e5SSatish Balay #undef __FUNC__
3255615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII"
3268f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
327416022c9SBarry Smith {
328416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
329496e697eSBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2;
33017ab2063SBarry Smith   FILE        *fd;
33117ab2063SBarry Smith   char        *outputname;
33217ab2063SBarry Smith 
3333a40ed3dSBarry Smith   PetscFunctionBegin;
33490ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
335416022c9SBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
33690ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
337a93ec695SBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
3383a40ed3dSBarry Smith     PetscFunctionReturn(0);
3393a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) {
340496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr);
341496e697eSBarry Smith     ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr);
342496e697eSBarry Smith     if (flg1 || flg2) fprintf(fd,"  not using I-node routines\n");
34395e01e2fSLois Curfman McInnes     else     fprintf(fd,"  using I-node routines: found %d nodes, limit used is %d\n",
34495e01e2fSLois Curfman McInnes         a->inode.node_count,a->inode.limit);
3453a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
346d00d2cf4SBarry Smith     int nofinalvalue = 0;
347d00d2cf4SBarry Smith     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) {
348d00d2cf4SBarry Smith       nofinalvalue = 1;
349d00d2cf4SBarry Smith     }
350416022c9SBarry Smith     fprintf(fd,"%% Size = %d %d \n",m,a->n);
3514e220ebcSLois Curfman McInnes     fprintf(fd,"%% Nonzeros = %d \n",a->nz);
352d00d2cf4SBarry Smith     fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue);
35317ab2063SBarry Smith     fprintf(fd,"zzz = [\n");
35417ab2063SBarry Smith 
35517ab2063SBarry Smith     for (i=0; i<m; i++) {
356416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
3573a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
358d64ed03dSBarry Smith         fprintf(fd,"%d %d  %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]),imag(a->a[j]));
35917ab2063SBarry Smith #else
3607a743949SBarry Smith         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j]+!shift, a->a[j]);
36117ab2063SBarry Smith #endif
36217ab2063SBarry Smith       }
36317ab2063SBarry Smith     }
364d00d2cf4SBarry Smith     if (nofinalvalue) {
365d00d2cf4SBarry Smith       fprintf(fd,"%d %d  %18.16e\n", m, a->n, 0.0);
366d00d2cf4SBarry Smith     }
36717ab2063SBarry Smith     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
3683a40ed3dSBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
36944cd7ae7SLois Curfman McInnes     for ( i=0; i<m; i++ ) {
37044cd7ae7SLois Curfman McInnes       fprintf(fd,"row %d:",i);
37144cd7ae7SLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
3723a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
373766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0 && real(a->a[j]) != 0.0)
37444cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
375766eeae4SLois Curfman McInnes         else if (imag(a->a[j]) < 0.0 && real(a->a[j]) != 0.0)
376766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
37744cd7ae7SLois Curfman McInnes         else if (real(a->a[j]) != 0.0)
37844cd7ae7SLois Curfman McInnes           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
37944cd7ae7SLois Curfman McInnes #else
38044cd7ae7SLois Curfman McInnes         if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
38144cd7ae7SLois Curfman McInnes #endif
38244cd7ae7SLois Curfman McInnes       }
38344cd7ae7SLois Curfman McInnes       fprintf(fd,"\n");
38444cd7ae7SLois Curfman McInnes     }
38544cd7ae7SLois Curfman McInnes   }
386496be53dSLois Curfman McInnes   else if (format == VIEWER_FORMAT_ASCII_SYMMODU) {
387496be53dSLois Curfman McInnes     int nzd=0, fshift=1, *sptr;
3882e44a96cSLois Curfman McInnes     sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr);
389496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
390496be53dSLois Curfman McInnes       sptr[i] = nzd+1;
391496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
392496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
3933a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
394496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) nzd++;
395496be53dSLois Curfman McInnes #else
396496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) nzd++;
397496be53dSLois Curfman McInnes #endif
398496be53dSLois Curfman McInnes         }
399496be53dSLois Curfman McInnes       }
400496be53dSLois Curfman McInnes     }
4012e44a96cSLois Curfman McInnes     sptr[m] = nzd+1;
402496be53dSLois Curfman McInnes     fprintf(fd," %d %d\n\n",m,nzd);
4032e44a96cSLois Curfman McInnes     for ( i=0; i<m+1; i+=6 ) {
4042e44a96cSLois 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]);
4052e44a96cSLois 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]);
4062e44a96cSLois 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]);
4072e44a96cSLois Curfman McInnes       else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]);
4082e44a96cSLois Curfman McInnes       else if (i<m)   fprintf(fd," %d %d\n",sptr[i],sptr[i+1]);
4097272d637SLois Curfman McInnes       else            fprintf(fd," %d\n",sptr[i]);
410496be53dSLois Curfman McInnes     }
411496be53dSLois Curfman McInnes     fprintf(fd,"\n");
412496be53dSLois Curfman McInnes     PetscFree(sptr);
413496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
414496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
415496be53dSLois Curfman McInnes         if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift);
416496be53dSLois Curfman McInnes       }
417496be53dSLois Curfman McInnes       fprintf(fd,"\n");
418496be53dSLois Curfman McInnes     }
419496be53dSLois Curfman McInnes     fprintf(fd,"\n");
420496be53dSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
421496be53dSLois Curfman McInnes       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
422496be53dSLois Curfman McInnes         if (a->j[j] >= i) {
4233a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
424496be53dSLois Curfman McInnes           if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0)
425496be53dSLois Curfman McInnes             fprintf(fd," %18.16e %18.16e ",real(a->a[j]),imag(a->a[j]));
426496be53dSLois Curfman McInnes #else
427496be53dSLois Curfman McInnes           if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]);
428496be53dSLois Curfman McInnes #endif
429496be53dSLois Curfman McInnes         }
430496be53dSLois Curfman McInnes       }
431496be53dSLois Curfman McInnes       fprintf(fd,"\n");
432496be53dSLois Curfman McInnes     }
4333a40ed3dSBarry Smith   } else {
43417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
43517ab2063SBarry Smith       fprintf(fd,"row %d:",i);
436416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
4373a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
438766eeae4SLois Curfman McInnes         if (imag(a->a[j]) > 0.0) {
439416022c9SBarry Smith           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
440766eeae4SLois Curfman McInnes         } else if (imag(a->a[j]) < 0.0) {
441766eeae4SLois Curfman McInnes           fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j]));
4423a40ed3dSBarry Smith         } else {
443416022c9SBarry Smith           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
44417ab2063SBarry Smith         }
44517ab2063SBarry Smith #else
446416022c9SBarry Smith         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
44717ab2063SBarry Smith #endif
44817ab2063SBarry Smith       }
44917ab2063SBarry Smith       fprintf(fd,"\n");
45017ab2063SBarry Smith     }
45117ab2063SBarry Smith   }
45217ab2063SBarry Smith   fflush(fd);
4533a40ed3dSBarry Smith   PetscFunctionReturn(0);
454416022c9SBarry Smith }
455416022c9SBarry Smith 
4565615d1e5SSatish Balay #undef __FUNC__
4575615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Draw"
4588f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
459416022c9SBarry Smith {
460416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
461cddf8d76SBarry Smith   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
4620513a670SBarry Smith   int         format;
46394a9d846SBarry Smith   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r,maxv = 0.0;
464bcd2baecSBarry Smith   Draw        draw;
465cddf8d76SBarry Smith   DrawButton  button;
46619bcc07fSBarry Smith   PetscTruth  isnull;
467cddf8d76SBarry Smith 
4683a40ed3dSBarry Smith   PetscFunctionBegin;
4690513a670SBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr);
47055843e3eSBarry Smith   ierr = DrawSynchronizedClear(draw); CHKERRQ(ierr);
4710513a670SBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
4723a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
47319bcc07fSBarry Smith 
474416022c9SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
475416022c9SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
476416022c9SBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
477416022c9SBarry Smith   /* loop over matrix elements drawing boxes */
4780513a670SBarry Smith 
4790513a670SBarry Smith   if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
4800513a670SBarry Smith     /* Blue for negative, Cyan for zero and  Red for positive */
481cddf8d76SBarry Smith     color = DRAW_BLUE;
482416022c9SBarry Smith     for ( i=0; i<m; i++ ) {
483cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
484416022c9SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
485cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
4863a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
487cddf8d76SBarry Smith         if (real(a->a[j]) >=  0.) continue;
488cddf8d76SBarry Smith #else
489cddf8d76SBarry Smith         if (a->a[j] >=  0.) continue;
490cddf8d76SBarry Smith #endif
491cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
492cddf8d76SBarry Smith       }
493cddf8d76SBarry Smith     }
494cddf8d76SBarry Smith     color = DRAW_CYAN;
495cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
496cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
497cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
498cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
499cddf8d76SBarry Smith         if (a->a[j] !=  0.) continue;
500cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
501cddf8d76SBarry Smith       }
502cddf8d76SBarry Smith     }
503cddf8d76SBarry Smith     color = DRAW_RED;
504cddf8d76SBarry Smith     for ( i=0; i<m; i++ ) {
505cddf8d76SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
506cddf8d76SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
507cddf8d76SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5083a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
509cddf8d76SBarry Smith         if (real(a->a[j]) <=  0.) continue;
510cddf8d76SBarry Smith #else
511cddf8d76SBarry Smith         if (a->a[j] <=  0.) continue;
512cddf8d76SBarry Smith #endif
513cddf8d76SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
514416022c9SBarry Smith       }
515416022c9SBarry Smith     }
5160513a670SBarry Smith   } else {
5170513a670SBarry Smith     /* use contour shading to indicate magnitude of values */
5180513a670SBarry Smith     /* first determine max of all nonzero values */
5190513a670SBarry Smith     int    nz = a->nz,count;
5200513a670SBarry Smith     Draw   popup;
5210513a670SBarry Smith 
5220513a670SBarry Smith     for ( i=0; i<nz; i++ ) {
5230513a670SBarry Smith       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
5240513a670SBarry Smith     }
525522c5e43SBarry Smith     ierr = DrawGetPopup(draw,&popup); CHKERRQ(ierr);
5260513a670SBarry Smith     ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr);
5270513a670SBarry Smith     count = 0;
5280513a670SBarry Smith     for ( i=0; i<m; i++ ) {
5290513a670SBarry Smith       y_l = m - i - 1.0; y_r = y_l + 1.0;
5300513a670SBarry Smith       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5310513a670SBarry Smith         x_l = a->j[j] + shift; x_r = x_l + 1.0;
5320513a670SBarry Smith         color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
5330513a670SBarry Smith         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
5340513a670SBarry Smith         count++;
5350513a670SBarry Smith       }
5360513a670SBarry Smith     }
5370513a670SBarry Smith   }
53855843e3eSBarry Smith   DrawSynchronizedFlush(draw);
539cddf8d76SBarry Smith   DrawGetPause(draw,&pause);
5403a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
541cddf8d76SBarry Smith 
542cddf8d76SBarry Smith   /* allow the matrix to zoom or shrink */
5436945ee14SBarry Smith   ierr = DrawCheckResizedWindow(draw);
54455843e3eSBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
545cddf8d76SBarry Smith   while (button != BUTTON_RIGHT) {
54655843e3eSBarry Smith     DrawSynchronizedClear(draw);
547cddf8d76SBarry Smith     if (button == BUTTON_LEFT) scale = .5;
548cddf8d76SBarry Smith     else if (button == BUTTON_CENTER) scale = 2.;
549cddf8d76SBarry Smith     xl = scale*(xl + w - xc) + xc - w*scale;
550cddf8d76SBarry Smith     xr = scale*(xr - w - xc) + xc + w*scale;
551cddf8d76SBarry Smith     yl = scale*(yl + h - yc) + yc - h*scale;
552cddf8d76SBarry Smith     yr = scale*(yr - h - yc) + yc + h*scale;
553cddf8d76SBarry Smith     w *= scale; h *= scale;
554cddf8d76SBarry Smith     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5550513a670SBarry Smith     if (format != VIEWER_FORMAT_DRAW_CONTOUR) {
5560513a670SBarry Smith       /* Blue for negative, Cyan for zero and  Red for positive */
557cddf8d76SBarry Smith       color = DRAW_BLUE;
558cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
559cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
560cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
561cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
5623a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
563cddf8d76SBarry Smith           if (real(a->a[j]) >=  0.) continue;
564cddf8d76SBarry Smith #else
565cddf8d76SBarry Smith           if (a->a[j] >=  0.) continue;
566cddf8d76SBarry Smith #endif
567cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
568cddf8d76SBarry Smith         }
569cddf8d76SBarry Smith       }
570cddf8d76SBarry Smith       color = DRAW_CYAN;
571cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
572cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
573cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
574cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
575cddf8d76SBarry Smith           if (a->a[j] !=  0.) continue;
576cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
577cddf8d76SBarry Smith         }
578cddf8d76SBarry Smith       }
579cddf8d76SBarry Smith       color = DRAW_RED;
580cddf8d76SBarry Smith       for ( i=0; i<m; i++ ) {
581cddf8d76SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
582cddf8d76SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
583cddf8d76SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
5843a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
585cddf8d76SBarry Smith           if (real(a->a[j]) <=  0.) continue;
586cddf8d76SBarry Smith #else
587cddf8d76SBarry Smith           if (a->a[j] <=  0.) continue;
588cddf8d76SBarry Smith #endif
589cddf8d76SBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
590cddf8d76SBarry Smith         }
591cddf8d76SBarry Smith       }
5920513a670SBarry Smith     } else {
5930513a670SBarry Smith       /* use contour shading to indicate magnitude of values */
5940513a670SBarry Smith       int count = 0;
5950513a670SBarry Smith       for ( i=0; i<m; i++ ) {
5960513a670SBarry Smith         y_l = m - i - 1.0; y_r = y_l + 1.0;
5970513a670SBarry Smith         for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
5980513a670SBarry Smith           x_l = a->j[j] + shift; x_r = x_l + 1.0;
5990513a670SBarry Smith           color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv);
6003a40ed3dSBarry Smith           DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); CHKERRQ(ierr);
6010513a670SBarry Smith           count++;
6020513a670SBarry Smith         }
6030513a670SBarry Smith       }
6040513a670SBarry Smith     }
6050513a670SBarry Smith 
6063a40ed3dSBarry Smith     ierr = DrawCheckResizedWindow(draw); CHKERRQ(ierr);
60755843e3eSBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);  CHKERRQ(ierr);
608cddf8d76SBarry Smith   }
6093a40ed3dSBarry Smith   PetscFunctionReturn(0);
610416022c9SBarry Smith }
611416022c9SBarry Smith 
6125615d1e5SSatish Balay #undef __FUNC__
613d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ"
614*e1311b90SBarry Smith int MatView_SeqAIJ(Mat A,Viewer viewer)
615416022c9SBarry Smith {
616416022c9SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
617bcd2baecSBarry Smith   ViewerType  vtype;
618bcd2baecSBarry Smith   int         ierr;
619416022c9SBarry Smith 
6203a40ed3dSBarry Smith   PetscFunctionBegin;
621bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
622bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
6233a40ed3dSBarry Smith     ierr = ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);  CHKERRQ(ierr);
624416022c9SBarry Smith   }
625bcd2baecSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
6263a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_ASCII(A,viewer); CHKERRQ(ierr);
627416022c9SBarry Smith   }
628bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
6293a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_Binary(A,viewer); CHKERRQ(ierr);
630416022c9SBarry Smith   }
631bcd2baecSBarry Smith   else if (vtype == DRAW_VIEWER) {
6323a40ed3dSBarry Smith     ierr = MatView_SeqAIJ_Draw(A,viewer); CHKERRQ(ierr);
63317ab2063SBarry Smith   }
6343a40ed3dSBarry Smith   PetscFunctionReturn(0);
63517ab2063SBarry Smith }
63619bcc07fSBarry Smith 
637c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat);
6385615d1e5SSatish Balay #undef __FUNC__
6395615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ"
6408f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
64117ab2063SBarry Smith {
642416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
64341c01911SSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr;
64443ee02c3SBarry Smith   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0;
645416022c9SBarry Smith   Scalar     *aa = a->a, *ap;
64617ab2063SBarry Smith 
6473a40ed3dSBarry Smith   PetscFunctionBegin;
6483a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
64917ab2063SBarry Smith 
65043ee02c3SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
65117ab2063SBarry Smith   for ( i=1; i<m; i++ ) {
652416022c9SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
65317ab2063SBarry Smith     fshift += imax[i-1] - ailen[i-1];
65494a9d846SBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
65517ab2063SBarry Smith     if (fshift) {
656416022c9SBarry Smith       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
65717ab2063SBarry Smith       N = ailen[i];
65817ab2063SBarry Smith       for ( j=0; j<N; j++ ) {
65917ab2063SBarry Smith         ip[j-fshift] = ip[j];
66017ab2063SBarry Smith         ap[j-fshift] = ap[j];
66117ab2063SBarry Smith       }
66217ab2063SBarry Smith     }
66317ab2063SBarry Smith     ai[i] = ai[i-1] + ailen[i-1];
66417ab2063SBarry Smith   }
66517ab2063SBarry Smith   if (m) {
66617ab2063SBarry Smith     fshift += imax[m-1] - ailen[m-1];
66717ab2063SBarry Smith     ai[m] = ai[m-1] + ailen[m-1];
66817ab2063SBarry Smith   }
66917ab2063SBarry Smith   /* reset ilen and imax for each row */
67017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
67117ab2063SBarry Smith     ailen[i] = imax[i] = ai[i+1] - ai[i];
67217ab2063SBarry Smith   }
673416022c9SBarry Smith   a->nz = ai[m] + shift;
67417ab2063SBarry Smith 
67517ab2063SBarry Smith   /* diagonals may have moved, so kill the diagonal pointers */
676416022c9SBarry Smith   if (fshift && a->diag) {
6770452661fSBarry Smith     PetscFree(a->diag);
678416022c9SBarry Smith     PLogObjectMemory(A,-(m+1)*sizeof(int));
679416022c9SBarry Smith     a->diag = 0;
68017ab2063SBarry Smith   }
6814e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n",
6824e220ebcSLois Curfman McInnes            m,a->n,fshift,a->nz);
6834e220ebcSLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n",
684b810aeb4SBarry Smith            a->reallocs);
68594a9d846SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax);
686dd5f02e7SSatish Balay   a->reallocs          = 0;
6874e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift;
6884e220ebcSLois Curfman McInnes 
68976dd722bSSatish Balay   /* check out for identical nodes. If found, use inode functions */
69041c01911SSatish Balay   ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr);
6913a40ed3dSBarry Smith   PetscFunctionReturn(0);
69217ab2063SBarry Smith }
69317ab2063SBarry Smith 
6945615d1e5SSatish Balay #undef __FUNC__
6955615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ"
6968f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A)
69717ab2063SBarry Smith {
698416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
6993a40ed3dSBarry Smith 
7003a40ed3dSBarry Smith   PetscFunctionBegin;
701cddf8d76SBarry Smith   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
7023a40ed3dSBarry Smith   PetscFunctionReturn(0);
70317ab2063SBarry Smith }
704416022c9SBarry Smith 
7055615d1e5SSatish Balay #undef __FUNC__
7065615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ"
707*e1311b90SBarry Smith int MatDestroy_SeqAIJ(Mat A)
70817ab2063SBarry Smith {
709416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
71082bf6240SBarry Smith   int        ierr;
711d5d45c9bSBarry Smith 
7123a40ed3dSBarry Smith   PetscFunctionBegin;
7133a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
714*e1311b90SBarry Smith   PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
71517ab2063SBarry Smith #endif
7160452661fSBarry Smith   PetscFree(a->a);
7170452661fSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
7180452661fSBarry Smith   if (a->diag) PetscFree(a->diag);
7190452661fSBarry Smith   if (a->ilen) PetscFree(a->ilen);
7200452661fSBarry Smith   if (a->imax) PetscFree(a->imax);
7210452661fSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
72276dd722bSSatish Balay   if (a->inode.size) PetscFree(a->inode.size);
72382bf6240SBarry Smith   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
7240452661fSBarry Smith   PetscFree(a);
725eed86810SBarry Smith 
726f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
727f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
7283a40ed3dSBarry Smith   PetscFunctionReturn(0);
72917ab2063SBarry Smith }
73017ab2063SBarry Smith 
7315615d1e5SSatish Balay #undef __FUNC__
7325615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ"
7338f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A)
73417ab2063SBarry Smith {
7353a40ed3dSBarry Smith   PetscFunctionBegin;
7363a40ed3dSBarry Smith   PetscFunctionReturn(0);
73717ab2063SBarry Smith }
73817ab2063SBarry Smith 
7395615d1e5SSatish Balay #undef __FUNC__
7405615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ"
7418f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op)
74217ab2063SBarry Smith {
743416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7443a40ed3dSBarry Smith 
7453a40ed3dSBarry Smith   PetscFunctionBegin;
7466d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
7476d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
7486d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
749219d9a1aSLois Curfman McInnes   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
7506d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
751c2653b3dSLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
75296854ed6SLois Curfman McInnes   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
7536d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
7546d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
755219d9a1aSLois Curfman McInnes            op == MAT_ROWS_UNSORTED ||
7566d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
7576d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
75890f02eecSBarry Smith            op == MAT_YES_NEW_DIAGONALS ||
759b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES||
760b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE)
761981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n");
7623a40ed3dSBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS) {
7633a40ed3dSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
7643a40ed3dSBarry Smith   } else if (op == MAT_INODE_LIMIT_1)            a->inode.limit  = 1;
7656d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
7666d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
7676d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
7686d4a8577SBarry Smith   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
7693a40ed3dSBarry Smith   else SETERRQ(PETSC_ERR_SUP,0,"unknown option");
7703a40ed3dSBarry Smith   PetscFunctionReturn(0);
77117ab2063SBarry Smith }
77217ab2063SBarry Smith 
7735615d1e5SSatish Balay #undef __FUNC__
7745615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ"
7758f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
77617ab2063SBarry Smith {
777416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
7783a40ed3dSBarry Smith   int        i,j, n,shift = a->indexshift,ierr;
77917ab2063SBarry Smith   Scalar     *x, zero = 0.0;
78017ab2063SBarry Smith 
7813a40ed3dSBarry Smith   PetscFunctionBegin;
7823a40ed3dSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
783*e1311b90SBarry Smith   ierr = VecGetArray(v,&x);;CHKERRQ(ierr);
784*e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr);
785a8c6a408SBarry Smith   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
786416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
787416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
788416022c9SBarry Smith       if (a->j[j]+shift == i) {
789416022c9SBarry Smith         x[i] = a->a[j];
79017ab2063SBarry Smith         break;
79117ab2063SBarry Smith       }
79217ab2063SBarry Smith     }
79317ab2063SBarry Smith   }
794*e1311b90SBarry Smith   ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr);
7953a40ed3dSBarry Smith   PetscFunctionReturn(0);
79617ab2063SBarry Smith }
79717ab2063SBarry Smith 
79817ab2063SBarry Smith /* -------------------------------------------------------*/
79917ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */
80017ab2063SBarry Smith /* -------------------------------------------------------*/
8015615d1e5SSatish Balay #undef __FUNC__
8025615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ"
80344cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
80417ab2063SBarry Smith {
805416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
80617ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
807*e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx, shift = a->indexshift;
80817ab2063SBarry Smith 
8093a40ed3dSBarry Smith   PetscFunctionBegin;
810*e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
811*e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
812cddf8d76SBarry Smith   PetscMemzero(y,a->n*sizeof(Scalar));
81317ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
81417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
815416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
816416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
817416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
81817ab2063SBarry Smith     alpha = x[i];
81917ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
82017ab2063SBarry Smith   }
821416022c9SBarry Smith   PLogFlops(2*a->nz - a->n);
822*e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
823*e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8243a40ed3dSBarry Smith   PetscFunctionReturn(0);
82517ab2063SBarry Smith }
826d5d45c9bSBarry Smith 
8275615d1e5SSatish Balay #undef __FUNC__
8285615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ"
82944cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
83017ab2063SBarry Smith {
831416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
83217ab2063SBarry Smith   Scalar     *x, *y, *v, alpha;
833*e1311b90SBarry Smith   int        ierr,m = a->m, n, i, *idx,shift = a->indexshift;
83417ab2063SBarry Smith 
8353a40ed3dSBarry Smith   PetscFunctionBegin;
836*e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
837*e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
83817ab2063SBarry Smith   if (zz != yy) VecCopy(zz,yy);
83917ab2063SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
84017ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
841416022c9SBarry Smith     idx   = a->j + a->i[i] + shift;
842416022c9SBarry Smith     v     = a->a + a->i[i] + shift;
843416022c9SBarry Smith     n     = a->i[i+1] - a->i[i];
84417ab2063SBarry Smith     alpha = x[i];
84517ab2063SBarry Smith     while (n-->0) {y[*idx++] += alpha * *v++;}
84617ab2063SBarry Smith   }
84790f02eecSBarry Smith   PLogFlops(2*a->nz);
848*e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
849*e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8503a40ed3dSBarry Smith   PetscFunctionReturn(0);
85117ab2063SBarry Smith }
85217ab2063SBarry Smith 
8535615d1e5SSatish Balay #undef __FUNC__
8545615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ"
85544cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
85617ab2063SBarry Smith {
857416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
85817ab2063SBarry Smith   Scalar     *x, *y, *v, sum;
859*e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
860e36a17ebSSatish Balay #if !defined(USE_FORTRAN_KERNELS)
861e36a17ebSSatish Balay   int        n, i, jrow,j;
862e36a17ebSSatish Balay #endif
86317ab2063SBarry Smith 
8643a40ed3dSBarry Smith   PetscFunctionBegin;
865*e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
866*e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
86717ab2063SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
868416022c9SBarry Smith   idx  = a->j;
869d64ed03dSBarry Smith   v    = a->a;
870416022c9SBarry Smith   ii   = a->i;
8718d195f9aSBarry Smith #if defined(USE_FORTRAN_KERNELS)
87229b5ca8aSSatish Balay   fortranmultaij_(&m,x,ii,idx+shift,v+shift,y);
8738d195f9aSBarry Smith #else
874d64ed03dSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
875d64ed03dSBarry Smith   idx  += shift;
87617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
8779ea0dfa2SSatish Balay     jrow = ii[i];
8789ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
87917ab2063SBarry Smith     sum  = 0.0;
8809ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
8819ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
8829ea0dfa2SSatish Balay      }
88317ab2063SBarry Smith     y[i] = sum;
88417ab2063SBarry Smith   }
8858d195f9aSBarry Smith #endif
886416022c9SBarry Smith   PLogFlops(2*a->nz - m);
887*e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
888*e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
8893a40ed3dSBarry Smith   PetscFunctionReturn(0);
89017ab2063SBarry Smith }
89117ab2063SBarry Smith 
8925615d1e5SSatish Balay #undef __FUNC__
8935615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ"
89444cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
89517ab2063SBarry Smith {
896416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
89717ab2063SBarry Smith   Scalar     *x, *y, *z, *v, sum;
898*e1311b90SBarry Smith   int        ierr,m = a->m, *idx, shift = a->indexshift,*ii;
899e36a17ebSSatish Balay #if !defined(USE_FORTRAN_KERNELS)
900e36a17ebSSatish Balay   int        n,i,jrow,j;
901e36a17ebSSatish Balay #endif
9029ea0dfa2SSatish Balay 
9033a40ed3dSBarry Smith   PetscFunctionBegin;
904*e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
905*e1311b90SBarry Smith   ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
906*e1311b90SBarry Smith   ierr = VecGetArray(zz,&z); CHKERRQ(ierr);
90717ab2063SBarry Smith   x    = x + shift; /* shift for Fortran start by 1 indexing */
908cddf8d76SBarry Smith   idx  = a->j;
909d64ed03dSBarry Smith   v    = a->a;
910cddf8d76SBarry Smith   ii   = a->i;
91102ab625aSSatish Balay #if defined(USE_FORTRAN_KERNELS)
91229b5ca8aSSatish Balay   fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
91302ab625aSSatish Balay #else
914d64ed03dSBarry Smith   v   += shift; /* shift for Fortran start by 1 indexing */
915d64ed03dSBarry Smith   idx += shift;
91617ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
9179ea0dfa2SSatish Balay     jrow = ii[i];
9189ea0dfa2SSatish Balay     n    = ii[i+1] - jrow;
91917ab2063SBarry Smith     sum  = y[i];
9209ea0dfa2SSatish Balay     for ( j=0; j<n; j++) {
9219ea0dfa2SSatish Balay       sum += v[jrow]*x[idx[jrow]]; jrow++;
9229ea0dfa2SSatish Balay      }
92317ab2063SBarry Smith     z[i] = sum;
92417ab2063SBarry Smith   }
92502ab625aSSatish Balay #endif
926416022c9SBarry Smith   PLogFlops(2*a->nz);
927*e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr);
928*e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr);
929*e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr);
9303a40ed3dSBarry Smith   PetscFunctionReturn(0);
93117ab2063SBarry Smith }
93217ab2063SBarry Smith 
93317ab2063SBarry Smith /*
93417ab2063SBarry Smith      Adds diagonal pointers to sparse matrix structure.
93517ab2063SBarry Smith */
93617ab2063SBarry Smith 
9375615d1e5SSatish Balay #undef __FUNC__
9385615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ"
939416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A)
94017ab2063SBarry Smith {
941416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
942416022c9SBarry Smith   int        i,j, *diag, m = a->m,shift = a->indexshift;
94317ab2063SBarry Smith 
9443a40ed3dSBarry Smith   PetscFunctionBegin;
9450452661fSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
946416022c9SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
947416022c9SBarry Smith   for ( i=0; i<a->m; i++ ) {
948416022c9SBarry Smith     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
949416022c9SBarry Smith       if (a->j[j]+shift == i) {
95017ab2063SBarry Smith         diag[i] = j - shift;
95117ab2063SBarry Smith         break;
95217ab2063SBarry Smith       }
95317ab2063SBarry Smith     }
95417ab2063SBarry Smith   }
955416022c9SBarry Smith   a->diag = diag;
9563a40ed3dSBarry Smith   PetscFunctionReturn(0);
95717ab2063SBarry Smith }
95817ab2063SBarry Smith 
9595615d1e5SSatish Balay #undef __FUNC__
9605615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ"
96144cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
96217ab2063SBarry Smith                            double fshift,int its,Vec xx)
96317ab2063SBarry Smith {
964416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
965416022c9SBarry Smith   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
966d5d45c9bSBarry Smith   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
96717ab2063SBarry Smith 
9683a40ed3dSBarry Smith   PetscFunctionBegin;
969*e1311b90SBarry Smith   ierr = VecGetArray(xx,&x); CHKERRQ(ierr);
970*e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
971d00d2cf4SBarry Smith   if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);}
972416022c9SBarry Smith   diag = a->diag;
973416022c9SBarry Smith   xs   = x + shift; /* shifted by one for index start of a or a->j*/
97417ab2063SBarry Smith   if (flag == SOR_APPLY_UPPER) {
97517ab2063SBarry Smith    /* apply ( U + D/omega) to the vector */
97617ab2063SBarry Smith     bs = b + shift;
97717ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
978416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
979416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
980416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
981416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
98217ab2063SBarry Smith         sum  = b[i]*d/omega;
98317ab2063SBarry Smith         SPARSEDENSEDOT(sum,bs,v,idx,n);
98417ab2063SBarry Smith         x[i] = sum;
98517ab2063SBarry Smith     }
9863a40ed3dSBarry Smith     PetscFunctionReturn(0);
98717ab2063SBarry Smith   }
98817ab2063SBarry Smith   if (flag == SOR_APPLY_LOWER) {
989a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done");
9903a40ed3dSBarry Smith   } else if (flag & SOR_EISENSTAT) {
99117ab2063SBarry Smith     /* Let  A = L + U + D; where L is lower trianglar,
99217ab2063SBarry Smith     U is upper triangular, E is diagonal; This routine applies
99317ab2063SBarry Smith 
99417ab2063SBarry Smith             (L + E)^{-1} A (U + E)^{-1}
99517ab2063SBarry Smith 
99617ab2063SBarry Smith     to a vector efficiently using Eisenstat's trick. This is for
99717ab2063SBarry Smith     the case of SSOR preconditioner, so E is D/omega where omega
99817ab2063SBarry Smith     is the relaxation factor.
99917ab2063SBarry Smith     */
10000452661fSBarry Smith     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
100117ab2063SBarry Smith     scale = (2.0/omega) - 1.0;
100217ab2063SBarry Smith 
100317ab2063SBarry Smith     /*  x = (E + U)^{-1} b */
100417ab2063SBarry Smith     for ( i=m-1; i>=0; i-- ) {
1005416022c9SBarry Smith       d    = fshift + a->a[diag[i] + shift];
1006416022c9SBarry Smith       n    = a->i[i+1] - diag[i] - 1;
1007416022c9SBarry Smith       idx  = a->j + diag[i] + (!shift);
1008416022c9SBarry Smith       v    = a->a + diag[i] + (!shift);
100917ab2063SBarry Smith       sum  = b[i];
101017ab2063SBarry Smith       SPARSEDENSEMDOT(sum,xs,v,idx,n);
101117ab2063SBarry Smith       x[i] = omega*(sum/d);
101217ab2063SBarry Smith     }
101317ab2063SBarry Smith 
101417ab2063SBarry Smith     /*  t = b - (2*E - D)x */
1015416022c9SBarry Smith     v = a->a;
101617ab2063SBarry Smith     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
101717ab2063SBarry Smith 
101817ab2063SBarry Smith     /*  t = (E + L)^{-1}t */
1019416022c9SBarry Smith     ts = t + shift; /* shifted by one for index start of a or a->j*/
1020416022c9SBarry Smith     diag = a->diag;
102117ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1022416022c9SBarry Smith       d    = fshift + a->a[diag[i]+shift];
1023416022c9SBarry Smith       n    = diag[i] - a->i[i];
1024416022c9SBarry Smith       idx  = a->j + a->i[i] + shift;
1025416022c9SBarry Smith       v    = a->a + a->i[i] + shift;
102617ab2063SBarry Smith       sum  = t[i];
102717ab2063SBarry Smith       SPARSEDENSEMDOT(sum,ts,v,idx,n);
102817ab2063SBarry Smith       t[i] = omega*(sum/d);
102917ab2063SBarry Smith     }
103017ab2063SBarry Smith 
103117ab2063SBarry Smith     /*  x = x + t */
103217ab2063SBarry Smith     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
10330452661fSBarry Smith     PetscFree(t);
10343a40ed3dSBarry Smith     PetscFunctionReturn(0);
103517ab2063SBarry Smith   }
103617ab2063SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
103717ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
103817ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1039416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1040416022c9SBarry Smith         n    = diag[i] - a->i[i];
1041416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1042416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
104317ab2063SBarry Smith         sum  = b[i];
104417ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
104517ab2063SBarry Smith         x[i] = omega*(sum/d);
104617ab2063SBarry Smith       }
104717ab2063SBarry Smith       xb = x;
10483a40ed3dSBarry Smith     } else xb = b;
104917ab2063SBarry Smith     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
105017ab2063SBarry Smith         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
105117ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1052416022c9SBarry Smith         x[i] *= a->a[diag[i]+shift];
105317ab2063SBarry Smith       }
105417ab2063SBarry Smith     }
105517ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
105617ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1057416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1058416022c9SBarry Smith         n    = a->i[i+1] - diag[i] - 1;
1059416022c9SBarry Smith         idx  = a->j + diag[i] + (!shift);
1060416022c9SBarry Smith         v    = a->a + diag[i] + (!shift);
106117ab2063SBarry Smith         sum  = xb[i];
106217ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
106317ab2063SBarry Smith         x[i] = omega*(sum/d);
106417ab2063SBarry Smith       }
106517ab2063SBarry Smith     }
106617ab2063SBarry Smith     its--;
106717ab2063SBarry Smith   }
106817ab2063SBarry Smith   while (its--) {
106917ab2063SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
107017ab2063SBarry Smith       for ( i=0; i<m; i++ ) {
1071416022c9SBarry Smith         d    = fshift + a->a[diag[i]+shift];
1072416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1073416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1074416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
107517ab2063SBarry Smith         sum  = b[i];
107617ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10777e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
107817ab2063SBarry Smith       }
107917ab2063SBarry Smith     }
108017ab2063SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
108117ab2063SBarry Smith       for ( i=m-1; i>=0; i-- ) {
1082416022c9SBarry Smith         d    = fshift + a->a[diag[i] + shift];
1083416022c9SBarry Smith         n    = a->i[i+1] - a->i[i];
1084416022c9SBarry Smith         idx  = a->j + a->i[i] + shift;
1085416022c9SBarry Smith         v    = a->a + a->i[i] + shift;
108617ab2063SBarry Smith         sum  = b[i];
108717ab2063SBarry Smith         SPARSEDENSEMDOT(sum,xs,v,idx,n);
10887e33a6baSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
108917ab2063SBarry Smith       }
109017ab2063SBarry Smith     }
109117ab2063SBarry Smith   }
10923a40ed3dSBarry Smith   PetscFunctionReturn(0);
109317ab2063SBarry Smith }
109417ab2063SBarry Smith 
10955615d1e5SSatish Balay #undef __FUNC__
10965615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ"
10978f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
109817ab2063SBarry Smith {
1099416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11004e220ebcSLois Curfman McInnes 
11013a40ed3dSBarry Smith   PetscFunctionBegin;
11024e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
11034e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
11044e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
11054e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
11064e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
11074e220ebcSLois Curfman McInnes   info->nz_allocated   = (double)a->maxnz;
11084e220ebcSLois Curfman McInnes   info->nz_used        = (double)a->nz;
11094e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(a->maxnz - a->nz);
11104e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
11114e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
11124e220ebcSLois Curfman McInnes   info->assemblies     = (double)A->num_ass;
11134e220ebcSLois Curfman McInnes   info->mallocs        = (double)a->reallocs;
11144e220ebcSLois Curfman McInnes   info->memory         = A->mem;
11154e220ebcSLois Curfman McInnes   if (A->factor) {
11164e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
11174e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
11184e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
11194e220ebcSLois Curfman McInnes   } else {
11204e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
11214e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
11224e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
11234e220ebcSLois Curfman McInnes   }
11243a40ed3dSBarry Smith   PetscFunctionReturn(0);
112517ab2063SBarry Smith }
112617ab2063SBarry Smith 
112717ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
112817ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
112917ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
113017ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
113117ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
113217ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
113317ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
113417ab2063SBarry Smith 
11355615d1e5SSatish Balay #undef __FUNC__
11365615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ"
11378f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
113817ab2063SBarry Smith {
1139416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1140416022c9SBarry Smith   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
114117ab2063SBarry Smith 
11423a40ed3dSBarry Smith   PetscFunctionBegin;
114377c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
114417ab2063SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
114517ab2063SBarry Smith   if (diag) {
114617ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1147a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1148416022c9SBarry Smith       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
1149416022c9SBarry Smith         a->ilen[rows[i]] = 1;
1150416022c9SBarry Smith         a->a[a->i[rows[i]]+shift] = *diag;
1151416022c9SBarry Smith         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
11523a40ed3dSBarry Smith       } else {
1153d64ed03dSBarry Smith         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr);
115417ab2063SBarry Smith       }
115517ab2063SBarry Smith     }
11563a40ed3dSBarry Smith   } else {
115717ab2063SBarry Smith     for ( i=0; i<N; i++ ) {
1158a8c6a408SBarry Smith       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
1159416022c9SBarry Smith       a->ilen[rows[i]] = 0;
116017ab2063SBarry Smith     }
116117ab2063SBarry Smith   }
116217ab2063SBarry Smith   ISRestoreIndices(is,&rows);
116343a90d84SBarry Smith   ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
11643a40ed3dSBarry Smith   PetscFunctionReturn(0);
116517ab2063SBarry Smith }
116617ab2063SBarry Smith 
11675615d1e5SSatish Balay #undef __FUNC__
11685615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ"
11698f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
117017ab2063SBarry Smith {
1171416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11723a40ed3dSBarry Smith 
11733a40ed3dSBarry Smith   PetscFunctionBegin;
11740752156aSBarry Smith   if (m) *m = a->m;
11750752156aSBarry Smith   if (n) *n = a->n;
11763a40ed3dSBarry Smith   PetscFunctionReturn(0);
117717ab2063SBarry Smith }
117817ab2063SBarry Smith 
11795615d1e5SSatish Balay #undef __FUNC__
11805615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ"
11818f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
118217ab2063SBarry Smith {
1183416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
11843a40ed3dSBarry Smith 
11853a40ed3dSBarry Smith   PetscFunctionBegin;
1186416022c9SBarry Smith   *m = 0; *n = a->m;
11873a40ed3dSBarry Smith   PetscFunctionReturn(0);
118817ab2063SBarry Smith }
1189026e39d0SSatish Balay 
11905615d1e5SSatish Balay #undef __FUNC__
11915615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ"
11924e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
119317ab2063SBarry Smith {
1194416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1195c456f294SBarry Smith   int        *itmp,i,shift = a->indexshift;
119617ab2063SBarry Smith 
11973a40ed3dSBarry Smith   PetscFunctionBegin;
1198a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
119917ab2063SBarry Smith 
1200416022c9SBarry Smith   *nz = a->i[row+1] - a->i[row];
1201416022c9SBarry Smith   if (v) *v = a->a + a->i[row] + shift;
120217ab2063SBarry Smith   if (idx) {
1203416022c9SBarry Smith     itmp = a->j + a->i[row] + shift;
12044e093b46SBarry Smith     if (*nz && shift) {
12050452661fSBarry Smith       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
120617ab2063SBarry Smith       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
12074e093b46SBarry Smith     } else if (*nz) {
12084e093b46SBarry Smith       *idx = itmp;
120917ab2063SBarry Smith     }
121017ab2063SBarry Smith     else *idx = 0;
121117ab2063SBarry Smith   }
12123a40ed3dSBarry Smith   PetscFunctionReturn(0);
121317ab2063SBarry Smith }
121417ab2063SBarry Smith 
12155615d1e5SSatish Balay #undef __FUNC__
12165615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ"
12174e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
121817ab2063SBarry Smith {
12194e093b46SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
12203a40ed3dSBarry Smith 
12213a40ed3dSBarry Smith   PetscFunctionBegin;
12224e093b46SBarry Smith   if (idx) {if (*idx && a->indexshift) PetscFree(*idx);}
12233a40ed3dSBarry Smith   PetscFunctionReturn(0);
122417ab2063SBarry Smith }
122517ab2063SBarry Smith 
12265615d1e5SSatish Balay #undef __FUNC__
12275615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ"
12288f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
122917ab2063SBarry Smith {
1230416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1231416022c9SBarry Smith   Scalar     *v = a->a;
123217ab2063SBarry Smith   double     sum = 0.0;
1233416022c9SBarry Smith   int        i, j,shift = a->indexshift;
123417ab2063SBarry Smith 
12353a40ed3dSBarry Smith   PetscFunctionBegin;
123617ab2063SBarry Smith   if (type == NORM_FROBENIUS) {
1237416022c9SBarry Smith     for (i=0; i<a->nz; i++ ) {
12383a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
123917ab2063SBarry Smith       sum += real(conj(*v)*(*v)); v++;
124017ab2063SBarry Smith #else
124117ab2063SBarry Smith       sum += (*v)*(*v); v++;
124217ab2063SBarry Smith #endif
124317ab2063SBarry Smith     }
124417ab2063SBarry Smith     *norm = sqrt(sum);
12453a40ed3dSBarry Smith   } else if (type == NORM_1) {
124617ab2063SBarry Smith     double *tmp;
1247416022c9SBarry Smith     int    *jj = a->j;
124866963ce1SSatish Balay     tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp);
1249cddf8d76SBarry Smith     PetscMemzero(tmp,a->n*sizeof(double));
125017ab2063SBarry Smith     *norm = 0.0;
1251416022c9SBarry Smith     for ( j=0; j<a->nz; j++ ) {
1252a2744918SBarry Smith         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
125317ab2063SBarry Smith     }
1254416022c9SBarry Smith     for ( j=0; j<a->n; j++ ) {
125517ab2063SBarry Smith       if (tmp[j] > *norm) *norm = tmp[j];
125617ab2063SBarry Smith     }
12570452661fSBarry Smith     PetscFree(tmp);
12583a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
125917ab2063SBarry Smith     *norm = 0.0;
1260416022c9SBarry Smith     for ( j=0; j<a->m; j++ ) {
1261416022c9SBarry Smith       v = a->a + a->i[j] + shift;
126217ab2063SBarry Smith       sum = 0.0;
1263416022c9SBarry Smith       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1264cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v++;
126517ab2063SBarry Smith       }
126617ab2063SBarry Smith       if (sum > *norm) *norm = sum;
126717ab2063SBarry Smith     }
12683a40ed3dSBarry Smith   } else {
1269a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"No support for two norm");
127017ab2063SBarry Smith   }
12713a40ed3dSBarry Smith   PetscFunctionReturn(0);
127217ab2063SBarry Smith }
127317ab2063SBarry Smith 
12745615d1e5SSatish Balay #undef __FUNC__
12755615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ"
12768f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B)
127717ab2063SBarry Smith {
1278416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
1279416022c9SBarry Smith   Mat        C;
1280416022c9SBarry Smith   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
1281416022c9SBarry Smith   int        shift = a->indexshift;
1282d5d45c9bSBarry Smith   Scalar     *array = a->a;
128317ab2063SBarry Smith 
12843a40ed3dSBarry Smith   PetscFunctionBegin;
1285a8c6a408SBarry Smith   if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place");
12860452661fSBarry Smith   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
1287cddf8d76SBarry Smith   PetscMemzero(col,(1+a->n)*sizeof(int));
128817ab2063SBarry Smith   if (shift) {
128917ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
129017ab2063SBarry Smith   }
129117ab2063SBarry Smith   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
1292416022c9SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
12930452661fSBarry Smith   PetscFree(col);
129417ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
129517ab2063SBarry Smith     len = ai[i+1]-ai[i];
1296416022c9SBarry Smith     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
129717ab2063SBarry Smith     array += len; aj += len;
129817ab2063SBarry Smith   }
129917ab2063SBarry Smith   if (shift) {
130017ab2063SBarry Smith     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
130117ab2063SBarry Smith   }
130217ab2063SBarry Smith 
13036d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
13046d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
130517ab2063SBarry Smith 
13063638b69dSLois Curfman McInnes   if (B != PETSC_NULL) {
1307416022c9SBarry Smith     *B = C;
130817ab2063SBarry Smith   } else {
1309f830108cSBarry Smith     PetscOps       *Abops;
1310f830108cSBarry Smith     struct _MatOps *Aops;
1311f830108cSBarry Smith 
1312416022c9SBarry Smith     /* This isn't really an in-place transpose */
13130452661fSBarry Smith     PetscFree(a->a);
13140452661fSBarry Smith     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
13150452661fSBarry Smith     if (a->diag) PetscFree(a->diag);
13160452661fSBarry Smith     if (a->ilen) PetscFree(a->ilen);
13170452661fSBarry Smith     if (a->imax) PetscFree(a->imax);
13180452661fSBarry Smith     if (a->solve_work) PetscFree(a->solve_work);
13191073c447SSatish Balay     if (a->inode.size) PetscFree(a->inode.size);
13200452661fSBarry Smith     PetscFree(a);
1321f830108cSBarry Smith 
1322f830108cSBarry Smith     /*
1323f830108cSBarry Smith         This is horrible, horrible code. We need to keep the
1324f830108cSBarry Smith       A pointers for the bops and ops but copy everything
1325f830108cSBarry Smith       else from C.
1326f830108cSBarry Smith     */
1327f830108cSBarry Smith     Abops = A->bops;
1328f830108cSBarry Smith     Aops  = A->ops;
1329f09e8eb9SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
1330f830108cSBarry Smith     A->bops = Abops;
1331f830108cSBarry Smith     A->ops  = Aops;
1332f830108cSBarry Smith 
13330452661fSBarry Smith     PetscHeaderDestroy(C);
133417ab2063SBarry Smith   }
13353a40ed3dSBarry Smith   PetscFunctionReturn(0);
133617ab2063SBarry Smith }
133717ab2063SBarry Smith 
13385615d1e5SSatish Balay #undef __FUNC__
13395615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ"
13408f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
134117ab2063SBarry Smith {
1342416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
134317ab2063SBarry Smith   Scalar     *l,*r,x,*v;
1344*e1311b90SBarry Smith   int        ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
134517ab2063SBarry Smith 
13463a40ed3dSBarry Smith   PetscFunctionBegin;
134717ab2063SBarry Smith   if (ll) {
13483ea7c6a1SSatish Balay     /* The local size is used so that VecMPI can be passed to this routine
13493ea7c6a1SSatish Balay        by MatDiagonalScale_MPIAIJ */
1350*e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr);
1351a8c6a408SBarry Smith     if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
1352*e1311b90SBarry Smith     ierr = VecGetArray(ll,&l); CHKERRQ(ierr);
1353416022c9SBarry Smith     v = a->a;
135417ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
135517ab2063SBarry Smith       x = l[i];
1356416022c9SBarry Smith       M = a->i[i+1] - a->i[i];
135717ab2063SBarry Smith       for ( j=0; j<M; j++ ) { (*v++) *= x;}
135817ab2063SBarry Smith     }
1359*e1311b90SBarry Smith     ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr);
136044cd7ae7SLois Curfman McInnes     PLogFlops(nz);
136117ab2063SBarry Smith   }
136217ab2063SBarry Smith   if (rr) {
1363*e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr);
1364a8c6a408SBarry Smith     if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
1365*e1311b90SBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1366416022c9SBarry Smith     v = a->a; jj = a->j;
136717ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
136817ab2063SBarry Smith       (*v++) *= r[*jj++ + shift];
136917ab2063SBarry Smith     }
1370*e1311b90SBarry Smith     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
137144cd7ae7SLois Curfman McInnes     PLogFlops(nz);
137217ab2063SBarry Smith   }
13733a40ed3dSBarry Smith   PetscFunctionReturn(0);
137417ab2063SBarry Smith }
137517ab2063SBarry Smith 
13765615d1e5SSatish Balay #undef __FUNC__
13775615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ"
13786a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B)
137917ab2063SBarry Smith {
1380db02288aSLois Curfman McInnes   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
1381d5db1b72SBarry Smith   int          *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
138299141d43SSatish Balay   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen;
1383a2744918SBarry Smith   register int sum,lensi;
138499141d43SSatish Balay   int          *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap;
138599141d43SSatish Balay   int          *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
138699141d43SSatish Balay   Scalar       *a_new,*mat_a;
1387416022c9SBarry Smith   Mat          C;
138817ab2063SBarry Smith 
13893a40ed3dSBarry Smith   PetscFunctionBegin;
1390d64ed03dSBarry Smith   ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr);
1391a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
1392d64ed03dSBarry Smith   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
1393a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
139499141d43SSatish Balay 
139517ab2063SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
139617ab2063SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
139717ab2063SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
139817ab2063SBarry Smith 
13997264ac53SSatish Balay   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */
140002834360SBarry Smith     /* special case of contiguous rows */
140157faeb66SBarry Smith     lens   = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens);
140202834360SBarry Smith     starts = lens + ncols;
140302834360SBarry Smith     /* loop over new rows determining lens and starting points */
140402834360SBarry Smith     for (i=0; i<nrows; i++) {
1405a2744918SBarry Smith       kstart  = ai[irow[i]]+shift;
1406a2744918SBarry Smith       kend    = kstart + ailen[irow[i]];
140702834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1408d8ced48eSBarry Smith         if (aj[k]+shift >= first) {
140902834360SBarry Smith           starts[i] = k;
141002834360SBarry Smith           break;
141102834360SBarry Smith 	}
141202834360SBarry Smith       }
1413a2744918SBarry Smith       sum = 0;
141402834360SBarry Smith       while (k < kend) {
1415d8ced48eSBarry Smith         if (aj[k++]+shift >= first+ncols) break;
1416a2744918SBarry Smith         sum++;
141702834360SBarry Smith       }
1418a2744918SBarry Smith       lens[i] = sum;
141902834360SBarry Smith     }
142002834360SBarry Smith     /* create submatrix */
1421cddf8d76SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
142208480c60SBarry Smith       int n_cols,n_rows;
142308480c60SBarry Smith       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1424a8c6a408SBarry Smith       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size");
1425d8ced48eSBarry Smith       ierr = MatZeroEntries(*B); CHKERRQ(ierr);
142608480c60SBarry Smith       C = *B;
14273a40ed3dSBarry Smith     } else {
142802834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
142908480c60SBarry Smith     }
1430db02288aSLois Curfman McInnes     c = (Mat_SeqAIJ*) C->data;
1431db02288aSLois Curfman McInnes 
143202834360SBarry Smith     /* loop over rows inserting into submatrix */
1433db02288aSLois Curfman McInnes     a_new    = c->a;
1434db02288aSLois Curfman McInnes     j_new    = c->j;
1435db02288aSLois Curfman McInnes     i_new    = c->i;
1436db02288aSLois Curfman McInnes     i_new[0] = -shift;
143702834360SBarry Smith     for (i=0; i<nrows; i++) {
1438a2744918SBarry Smith       ii    = starts[i];
1439a2744918SBarry Smith       lensi = lens[i];
1440a2744918SBarry Smith       for ( k=0; k<lensi; k++ ) {
1441a2744918SBarry Smith         *j_new++ = aj[ii+k] - first;
144202834360SBarry Smith       }
1443a2744918SBarry Smith       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1444a2744918SBarry Smith       a_new      += lensi;
1445a2744918SBarry Smith       i_new[i+1]  = i_new[i] + lensi;
1446a2744918SBarry Smith       c->ilen[i]  = lensi;
144702834360SBarry Smith     }
14480452661fSBarry Smith     PetscFree(lens);
14493a40ed3dSBarry Smith   } else {
145002834360SBarry Smith     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
14510452661fSBarry Smith     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
145202834360SBarry Smith     ssmap = smap + shift;
145399141d43SSatish Balay     lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
1454cddf8d76SBarry Smith     PetscMemzero(smap,oldcols*sizeof(int));
145517ab2063SBarry Smith     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
145602834360SBarry Smith     /* determine lens of each row */
145702834360SBarry Smith     for (i=0; i<nrows; i++) {
1458d8ced48eSBarry Smith       kstart  = ai[irow[i]]+shift;
145902834360SBarry Smith       kend    = kstart + a->ilen[irow[i]];
146002834360SBarry Smith       lens[i] = 0;
146102834360SBarry Smith       for ( k=kstart; k<kend; k++ ) {
1462d8ced48eSBarry Smith         if (ssmap[aj[k]]) {
146302834360SBarry Smith           lens[i]++;
146402834360SBarry Smith         }
146502834360SBarry Smith       }
146602834360SBarry Smith     }
146717ab2063SBarry Smith     /* Create and fill new matrix */
1468a2744918SBarry Smith     if (scall == MAT_REUSE_MATRIX) {
146999141d43SSatish Balay       c = (Mat_SeqAIJ *)((*B)->data);
147099141d43SSatish Balay 
1471a8c6a408SBarry Smith       if (c->m  != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
147299141d43SSatish Balay       if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) {
1473a8c6a408SBarry Smith         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
147499141d43SSatish Balay       }
147599141d43SSatish Balay       PetscMemzero(c->ilen,c->m*sizeof(int));
147608480c60SBarry Smith       C = *B;
14773a40ed3dSBarry Smith     } else {
147802834360SBarry Smith       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
147908480c60SBarry Smith     }
148099141d43SSatish Balay     c = (Mat_SeqAIJ *)(C->data);
148117ab2063SBarry Smith     for (i=0; i<nrows; i++) {
148299141d43SSatish Balay       row    = irow[i];
148399141d43SSatish Balay       kstart = ai[row]+shift;
148499141d43SSatish Balay       kend   = kstart + a->ilen[row];
148599141d43SSatish Balay       mat_i  = c->i[i]+shift;
148699141d43SSatish Balay       mat_j  = c->j + mat_i;
148799141d43SSatish Balay       mat_a  = c->a + mat_i;
148899141d43SSatish Balay       mat_ilen = c->ilen + i;
148917ab2063SBarry Smith       for ( k=kstart; k<kend; k++ ) {
149099141d43SSatish Balay         if ((tcol=ssmap[a->j[k]])) {
149199141d43SSatish Balay           *mat_j++ = tcol - (!shift);
149299141d43SSatish Balay           *mat_a++ = a->a[k];
149399141d43SSatish Balay           (*mat_ilen)++;
149499141d43SSatish Balay 
149517ab2063SBarry Smith         }
149617ab2063SBarry Smith       }
149717ab2063SBarry Smith     }
149802834360SBarry Smith     /* Free work space */
149902834360SBarry Smith     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
150099141d43SSatish Balay     PetscFree(smap); PetscFree(lens);
150102834360SBarry Smith   }
15026d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
15036d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
150417ab2063SBarry Smith 
150517ab2063SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1506416022c9SBarry Smith   *B = C;
15073a40ed3dSBarry Smith   PetscFunctionReturn(0);
150817ab2063SBarry Smith }
150917ab2063SBarry Smith 
1510a871dcd8SBarry Smith /*
151163b91edcSBarry Smith      note: This can only work for identity for row and col. It would
151263b91edcSBarry Smith    be good to check this and otherwise generate an error.
1513a871dcd8SBarry Smith */
15145615d1e5SSatish Balay #undef __FUNC__
15155615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ"
15168f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1517a871dcd8SBarry Smith {
151863b91edcSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
151908480c60SBarry Smith   int        ierr;
152063b91edcSBarry Smith   Mat        outA;
152163b91edcSBarry Smith 
15223a40ed3dSBarry Smith   PetscFunctionBegin;
1523a8c6a408SBarry Smith   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
1524a871dcd8SBarry Smith 
152563b91edcSBarry Smith   outA          = inA;
152663b91edcSBarry Smith   inA->factor   = FACTOR_LU;
152763b91edcSBarry Smith   a->row        = row;
152863b91edcSBarry Smith   a->col        = col;
152963b91edcSBarry Smith 
1530f0ec6fceSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1531f0ec6fceSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
15321e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1533f0ec6fceSSatish Balay 
153494a9d846SBarry Smith   if (!a->solve_work) { /* this matrix may have been factored before */
15350452661fSBarry Smith     a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
153694a9d846SBarry Smith   }
153763b91edcSBarry Smith 
153808480c60SBarry Smith   if (!a->diag) {
153908480c60SBarry Smith     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
154063b91edcSBarry Smith   }
154163b91edcSBarry Smith   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
15423a40ed3dSBarry Smith   PetscFunctionReturn(0);
1543a871dcd8SBarry Smith }
1544a871dcd8SBarry Smith 
1545f0b747eeSBarry Smith #include "pinclude/plapack.h"
15465615d1e5SSatish Balay #undef __FUNC__
15475615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ"
15488f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1549f0b747eeSBarry Smith {
1550f0b747eeSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1551f0b747eeSBarry Smith   int        one = 1;
15523a40ed3dSBarry Smith 
15533a40ed3dSBarry Smith   PetscFunctionBegin;
1554f0b747eeSBarry Smith   BLscal_( &a->nz, alpha, a->a, &one );
1555f0b747eeSBarry Smith   PLogFlops(a->nz);
15563a40ed3dSBarry Smith   PetscFunctionReturn(0);
1557f0b747eeSBarry Smith }
1558f0b747eeSBarry Smith 
15595615d1e5SSatish Balay #undef __FUNC__
15605615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ"
15616a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
1562cddf8d76SBarry Smith {
1563cddf8d76SBarry Smith   int ierr,i;
1564cddf8d76SBarry Smith 
15653a40ed3dSBarry Smith   PetscFunctionBegin;
1566cddf8d76SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
15670452661fSBarry Smith     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1568cddf8d76SBarry Smith   }
1569cddf8d76SBarry Smith 
1570cddf8d76SBarry Smith   for ( i=0; i<n; i++ ) {
15716a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1572cddf8d76SBarry Smith   }
15733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1574cddf8d76SBarry Smith }
1575cddf8d76SBarry Smith 
15765615d1e5SSatish Balay #undef __FUNC__
15775615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ"
15788f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs)
15795a838052SSatish Balay {
1580f830108cSBarry Smith   PetscFunctionBegin;
15815a838052SSatish Balay   *bs = 1;
15823a40ed3dSBarry Smith   PetscFunctionReturn(0);
15835a838052SSatish Balay }
15845a838052SSatish Balay 
15855615d1e5SSatish Balay #undef __FUNC__
15865615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ"
15878f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov)
15884dcbc457SBarry Smith {
1589e4d965acSSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
159006763907SSatish Balay   int        shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val;
15918a047759SSatish Balay   int        start, end, *ai, *aj;
1592bbd702dbSSatish Balay   BT         table;
1593bbd702dbSSatish Balay 
15943a40ed3dSBarry Smith   PetscFunctionBegin;
15958a047759SSatish Balay   shift = a->indexshift;
1596e4d965acSSatish Balay   m     = a->m;
1597e4d965acSSatish Balay   ai    = a->i;
15988a047759SSatish Balay   aj    = a->j+shift;
15998a047759SSatish Balay 
1600a8c6a408SBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used");
160106763907SSatish Balay 
160206763907SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
1603bbd702dbSSatish Balay   ierr  = BTCreate(m,table); CHKERRQ(ierr);
160406763907SSatish Balay 
1605e4d965acSSatish Balay   for ( i=0; i<is_max; i++ ) {
1606b97fc60eSLois Curfman McInnes     /* Initialize the two local arrays */
1607e4d965acSSatish Balay     isz  = 0;
1608bbd702dbSSatish Balay     BTMemzero(m,table);
1609e4d965acSSatish Balay 
1610e4d965acSSatish Balay     /* Extract the indices, assume there can be duplicate entries */
16114dcbc457SBarry Smith     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
161277c4ece6SBarry Smith     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
1613e4d965acSSatish Balay 
1614dd097bc3SLois Curfman McInnes     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1615e4d965acSSatish Balay     for ( j=0; j<n ; ++j){
1616bbd702dbSSatish Balay       if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];}
16174dcbc457SBarry Smith     }
161806763907SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
161906763907SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
1620e4d965acSSatish Balay 
162104a348a9SBarry Smith     k = 0;
162204a348a9SBarry Smith     for ( j=0; j<ov; j++){ /* for each overlap */
162304a348a9SBarry Smith       n = isz;
162406763907SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1625e4d965acSSatish Balay         row   = nidx[k];
1626e4d965acSSatish Balay         start = ai[row];
1627e4d965acSSatish Balay         end   = ai[row+1];
162804a348a9SBarry Smith         for ( l = start; l<end ; l++){
16298a047759SSatish Balay           val = aj[l] + shift;
16302a8f2162SSatish Balay           if (!BTLookupSet(table,val)) {nidx[isz++] = val;}
1631e4d965acSSatish Balay         }
1632e4d965acSSatish Balay       }
1633e4d965acSSatish Balay     }
1634029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr);
1635e4d965acSSatish Balay   }
1636bbd702dbSSatish Balay   BTDestroy(table);
163706763907SSatish Balay   PetscFree(nidx);
16383a40ed3dSBarry Smith   PetscFunctionReturn(0);
16394dcbc457SBarry Smith }
164017ab2063SBarry Smith 
16410513a670SBarry Smith /* -------------------------------------------------------------- */
16420513a670SBarry Smith #undef __FUNC__
16430513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ"
16440513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B)
16450513a670SBarry Smith {
16460513a670SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
16470513a670SBarry Smith   Scalar     *vwork;
16480513a670SBarry Smith   int        i, ierr, nz, m = a->m, n = a->n, *cwork;
16490513a670SBarry Smith   int        *row,*col,*cnew,j,*lens;
165056cd22aeSBarry Smith   IS         icolp,irowp;
16510513a670SBarry Smith 
16523a40ed3dSBarry Smith   PetscFunctionBegin;
165356cd22aeSBarry Smith   ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr);
165456cd22aeSBarry Smith   ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr);
165556cd22aeSBarry Smith   ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr);
165656cd22aeSBarry Smith   ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr);
16570513a670SBarry Smith 
16580513a670SBarry Smith   /* determine lengths of permuted rows */
16590513a670SBarry Smith   lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens);
16600513a670SBarry Smith   for (i=0; i<m; i++ ) {
16610513a670SBarry Smith     lens[row[i]] = a->i[i+1] - a->i[i];
16620513a670SBarry Smith   }
16630513a670SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr);
16640513a670SBarry Smith   PetscFree(lens);
16650513a670SBarry Smith 
16660513a670SBarry Smith   cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew);
16670513a670SBarry Smith   for (i=0; i<m; i++) {
16680513a670SBarry Smith     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16690513a670SBarry Smith     for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];}
16700513a670SBarry Smith     ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr);
16710513a670SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
16720513a670SBarry Smith   }
16730513a670SBarry Smith   PetscFree(cnew);
16740513a670SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
16750513a670SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
167656cd22aeSBarry Smith   ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr);
167756cd22aeSBarry Smith   ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr);
167856cd22aeSBarry Smith   ierr = ISDestroy(irowp); CHKERRQ(ierr);
167956cd22aeSBarry Smith   ierr = ISDestroy(icolp); CHKERRQ(ierr);
16803a40ed3dSBarry Smith   PetscFunctionReturn(0);
16810513a670SBarry Smith }
16820513a670SBarry Smith 
16835615d1e5SSatish Balay #undef __FUNC__
16845615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ"
1685682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A)
1686682d7d0cSBarry Smith {
1687682d7d0cSBarry Smith   static int called = 0;
1688682d7d0cSBarry Smith   MPI_Comm   comm = A->comm;
1689682d7d0cSBarry Smith 
16903a40ed3dSBarry Smith   PetscFunctionBegin;
16913a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
169276be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
169376be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
169476be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
169576be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodes\n");
169676be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1697682d7d0cSBarry Smith #if defined(HAVE_ESSL)
169876be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.\n");
1699682d7d0cSBarry Smith #endif
17003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1701682d7d0cSBarry Smith }
17028f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg);
1703a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1704a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *);
1705a93ec695SBarry Smith 
1706682d7d0cSBarry Smith /* -------------------------------------------------------------------*/
170717ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
170817ab2063SBarry Smith        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1709416022c9SBarry Smith        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1710416022c9SBarry Smith        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
171117ab2063SBarry Smith        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
171217ab2063SBarry Smith        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
171317ab2063SBarry Smith        MatLUFactor_SeqAIJ,0,
171417ab2063SBarry Smith        MatRelax_SeqAIJ,
171517ab2063SBarry Smith        MatTranspose_SeqAIJ,
17167264ac53SSatish Balay        MatGetInfo_SeqAIJ,MatEqual_SeqAIJ,
1717f0b747eeSBarry Smith        MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ,
171817ab2063SBarry Smith        0,MatAssemblyEnd_SeqAIJ,
171917ab2063SBarry Smith        MatCompress_SeqAIJ,
172017ab2063SBarry Smith        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
172117ab2063SBarry Smith        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
172217ab2063SBarry Smith        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
172317ab2063SBarry Smith        MatILUFactorSymbolic_SeqAIJ,0,
172494a9d846SBarry Smith        0,0,
17253d1612f7SBarry Smith        MatConvertSameType_SeqAIJ,0,0,
1726cddf8d76SBarry Smith        MatILUFactor_SeqAIJ,0,0,
17277eb43aa7SLois Curfman McInnes        MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ,
1728682d7d0cSBarry Smith        MatGetValues_SeqAIJ,0,
1729f0b747eeSBarry Smith        MatPrintHelp_SeqAIJ,
17305a838052SSatish Balay        MatScale_SeqAIJ,0,0,
17316945ee14SBarry Smith        MatILUDTFactor_SeqAIJ,
17326945ee14SBarry Smith        MatGetBlockSize_SeqAIJ,
17333b2fbd54SBarry Smith        MatGetRowIJ_SeqAIJ,
17343b2fbd54SBarry Smith        MatRestoreRowIJ_SeqAIJ,
17353b2fbd54SBarry Smith        MatGetColumnIJ_SeqAIJ,
1736a93ec695SBarry Smith        MatRestoreColumnIJ_SeqAIJ,
1737a93ec695SBarry Smith        MatFDColoringCreate_SeqAIJ,
17380513a670SBarry Smith        MatColoringPatch_SeqAIJ,
17390513a670SBarry Smith        0,
17400513a670SBarry Smith        MatPermute_SeqAIJ};
174117ab2063SBarry Smith 
174217ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat);
174317ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat);
174417ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat);
174517ab2063SBarry Smith 
17465615d1e5SSatish Balay #undef __FUNC__
17475615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ"
174817ab2063SBarry Smith /*@C
1749682d7d0cSBarry Smith    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
17500d15e28bSLois Curfman McInnes    (the default parallel PETSc format).  For good matrix assembly performance
17516e62573dSLois Curfman McInnes    the user should preallocate the matrix storage by setting the parameter nz
17522bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
17532bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
175417ab2063SBarry Smith 
175517ab2063SBarry Smith    Input Parameters:
1756029af93fSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
175717ab2063SBarry Smith .  m - number of rows
175817ab2063SBarry Smith .  n - number of columns
175917ab2063SBarry Smith .  nz - number of nonzeros per row (same for all rows)
17602bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of nonzeros in the various rows
17612bd5e0b2SLois Curfman McInnes          (possibly different for each row) or PETSC_NULL
176217ab2063SBarry Smith 
176317ab2063SBarry Smith    Output Parameter:
1764416022c9SBarry Smith .  A - the matrix
176517ab2063SBarry Smith 
176617ab2063SBarry Smith    Notes:
176717ab2063SBarry Smith    The AIJ format (also called the Yale sparse matrix format or
176817ab2063SBarry Smith    compressed row storage), is fully compatible with standard Fortran 77
17690002213bSLois Curfman McInnes    storage.  That is, the stored row and column indices can begin at
177044cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
177117ab2063SBarry Smith 
177217ab2063SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1773a40aa06bSLois Curfman McInnes    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
17743d323bbdSBarry Smith    allocation.  For large problems you MUST preallocate memory or you
17756da5968aSLois Curfman McInnes    will get TERRIBLE performance, see the users' manual chapter on matrices.
177617ab2063SBarry Smith 
1777682d7d0cSBarry Smith    By default, this format uses inodes (identical nodes) when possible, to
17784fca80b9SLois Curfman McInnes    improve numerical efficiency of matrix-vector products and solves. We
1779682d7d0cSBarry Smith    search for consecutive rows with the same nonzero structure, thereby
17806c7ebb05SLois Curfman McInnes    reusing matrix information to achieve increased efficiency.
17816c7ebb05SLois Curfman McInnes 
17826c7ebb05SLois Curfman McInnes    Options Database Keys:
17836c7ebb05SLois Curfman McInnes $    -mat_aij_no_inode  - Do not use inodes
17846e62573dSLois Curfman McInnes $    -mat_aij_inode_limit <limit> - Set inode limit.
17856e62573dSLois Curfman McInnes $        (max limit=5)
17866e62573dSLois Curfman McInnes $    -mat_aij_oneindex - Internally use indexing starting at 1
17876e62573dSLois Curfman McInnes $        rather than 0.  Note: When calling MatSetValues(),
17886e62573dSLois Curfman McInnes $        the user still MUST index entries starting at 0!
178917ab2063SBarry Smith 
179017ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
179117ab2063SBarry Smith @*/
1792416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
179317ab2063SBarry Smith {
1794416022c9SBarry Smith   Mat        B;
1795416022c9SBarry Smith   Mat_SeqAIJ *b;
17966945ee14SBarry Smith   int        i, len, ierr, flg,size;
17976945ee14SBarry Smith 
17983a40ed3dSBarry Smith   PetscFunctionBegin;
17996945ee14SBarry Smith   MPI_Comm_size(comm,&size);
1800a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1");
1801d5d45c9bSBarry Smith 
1802416022c9SBarry Smith   *A                  = 0;
1803f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView);
1804416022c9SBarry Smith   PLogObjectCreate(B);
18050452661fSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
180644cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqAIJ));
1807f830108cSBarry Smith   PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps));
1808*e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_SeqAIJ;
1809*e1311b90SBarry Smith   B->ops->view             = MatView_SeqAIJ;
1810416022c9SBarry Smith   B->factor           = 0;
1811416022c9SBarry Smith   B->lupivotthreshold = 1.0;
181290f02eecSBarry Smith   B->mapping          = 0;
1813*e1311b90SBarry Smith   ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg); CHKERRQ(ierr);
18147a743949SBarry Smith   b->ilu_preserve_row_sums = PETSC_FALSE;
1815*e1311b90SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*) &b->ilu_preserve_row_sums);CHKERRQ(ierr);
1816416022c9SBarry Smith   b->row              = 0;
1817416022c9SBarry Smith   b->col              = 0;
181882bf6240SBarry Smith   b->icol             = 0;
1819416022c9SBarry Smith   b->indexshift       = 0;
1820b810aeb4SBarry Smith   b->reallocs         = 0;
182169957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr);
182269957df2SSatish Balay   if (flg) b->indexshift = -1;
182317ab2063SBarry Smith 
182444cd7ae7SLois Curfman McInnes   b->m = m; B->m = m; B->M = m;
182544cd7ae7SLois Curfman McInnes   b->n = n; B->n = n; B->N = n;
18260452661fSBarry Smith   b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1827b4fd4287SBarry Smith   if (nnz == PETSC_NULL) {
18287b8455f0SLois Curfman McInnes     if (nz == PETSC_DEFAULT) nz = 10;
18297b8455f0SLois Curfman McInnes     else if (nz <= 0)        nz = 1;
1830416022c9SBarry Smith     for ( i=0; i<m; i++ ) b->imax[i] = nz;
183117ab2063SBarry Smith     nz = nz*m;
18323a40ed3dSBarry Smith   } else {
183317ab2063SBarry Smith     nz = 0;
1834416022c9SBarry Smith     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
183517ab2063SBarry Smith   }
183617ab2063SBarry Smith 
183717ab2063SBarry Smith   /* allocate the matrix space */
1838416022c9SBarry Smith   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
18390452661fSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1840416022c9SBarry Smith   b->j  = (int *) (b->a + nz);
1841cddf8d76SBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
1842416022c9SBarry Smith   b->i  = b->j + nz;
1843416022c9SBarry Smith   b->singlemalloc = 1;
184417ab2063SBarry Smith 
1845416022c9SBarry Smith   b->i[0] = -b->indexshift;
184617ab2063SBarry Smith   for (i=1; i<m+1; i++) {
1847416022c9SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
184817ab2063SBarry Smith   }
184917ab2063SBarry Smith 
1850416022c9SBarry Smith   /* b->ilen will count nonzeros in each row so far. */
18510452661fSBarry Smith   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1852f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1853416022c9SBarry Smith   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
185417ab2063SBarry Smith 
1855416022c9SBarry Smith   b->nz               = 0;
1856416022c9SBarry Smith   b->maxnz            = nz;
1857416022c9SBarry Smith   b->sorted           = 0;
1858416022c9SBarry Smith   b->roworiented      = 1;
1859416022c9SBarry Smith   b->nonew            = 0;
1860416022c9SBarry Smith   b->diag             = 0;
1861416022c9SBarry Smith   b->solve_work       = 0;
186271bd300dSLois Curfman McInnes   b->spptr            = 0;
1863754ec7b1SSatish Balay   b->inode.node_count = 0;
1864754ec7b1SSatish Balay   b->inode.size       = 0;
18656c7ebb05SLois Curfman McInnes   b->inode.limit      = 5;
18666c7ebb05SLois Curfman McInnes   b->inode.max_limit  = 5;
18674e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
186817ab2063SBarry Smith 
1869416022c9SBarry Smith   *A = B;
18704e220ebcSLois Curfman McInnes 
18714b14c69eSBarry Smith   /*  SuperLU is not currently supported through PETSc */
18724b14c69eSBarry Smith #if defined(HAVE_SUPERLU)
187369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr);
187469957df2SSatish Balay   if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); }
18754b14c69eSBarry Smith #endif
187669957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr);
187769957df2SSatish Balay   if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); }
187869957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr);
187969957df2SSatish Balay   if (flg) {
1880a8c6a408SBarry Smith     if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml");
1881416022c9SBarry Smith     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
188217ab2063SBarry Smith   }
188369957df2SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
188469957df2SSatish Balay   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
18853a40ed3dSBarry Smith   PetscFunctionReturn(0);
188617ab2063SBarry Smith }
188717ab2063SBarry Smith 
18885615d1e5SSatish Balay #undef __FUNC__
18895615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ"
18903d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues)
189117ab2063SBarry Smith {
1892416022c9SBarry Smith   Mat        C;
1893416022c9SBarry Smith   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
189408480c60SBarry Smith   int        i,len, m = a->m,shift = a->indexshift;
189517ab2063SBarry Smith 
18963a40ed3dSBarry Smith   PetscFunctionBegin;
18974043dd9cSLois Curfman McInnes   *B = 0;
1898f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView);
1899416022c9SBarry Smith   PLogObjectCreate(C);
19000452661fSBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1901f830108cSBarry Smith   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1902*e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqAIJ;
1903*e1311b90SBarry Smith   C->ops->view       = MatView_SeqAIJ;
1904416022c9SBarry Smith   C->factor     = A->factor;
1905416022c9SBarry Smith   c->row        = 0;
1906416022c9SBarry Smith   c->col        = 0;
190782bf6240SBarry Smith   c->icol       = 0;
1908416022c9SBarry Smith   c->indexshift = shift;
1909c456f294SBarry Smith   C->assembled  = PETSC_TRUE;
191017ab2063SBarry Smith 
191144cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
191244cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
191344cd7ae7SLois Curfman McInnes   C->M          = a->m;
191444cd7ae7SLois Curfman McInnes   C->N          = a->n;
191517ab2063SBarry Smith 
19160452661fSBarry Smith   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
19170452661fSBarry Smith   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
191817ab2063SBarry Smith   for ( i=0; i<m; i++ ) {
1919416022c9SBarry Smith     c->imax[i] = a->imax[i];
1920416022c9SBarry Smith     c->ilen[i] = a->ilen[i];
192117ab2063SBarry Smith   }
192217ab2063SBarry Smith 
192317ab2063SBarry Smith   /* allocate the matrix space */
1924416022c9SBarry Smith   c->singlemalloc = 1;
1925416022c9SBarry Smith   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
19260452661fSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1927416022c9SBarry Smith   c->j  = (int *) (c->a + a->i[m] + shift);
1928416022c9SBarry Smith   c->i  = c->j + a->i[m] + shift;
1929416022c9SBarry Smith   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
193017ab2063SBarry Smith   if (m > 0) {
1931416022c9SBarry Smith     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
193208480c60SBarry Smith     if (cpvalues == COPY_VALUES) {
1933416022c9SBarry Smith       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
193417ab2063SBarry Smith     }
193508480c60SBarry Smith   }
193617ab2063SBarry Smith 
1937f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
1938416022c9SBarry Smith   c->sorted      = a->sorted;
1939416022c9SBarry Smith   c->roworiented = a->roworiented;
1940416022c9SBarry Smith   c->nonew       = a->nonew;
19417a743949SBarry Smith   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
194217ab2063SBarry Smith 
1943416022c9SBarry Smith   if (a->diag) {
19440452661fSBarry Smith     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1945416022c9SBarry Smith     PLogObjectMemory(C,(m+1)*sizeof(int));
194617ab2063SBarry Smith     for ( i=0; i<m; i++ ) {
1947416022c9SBarry Smith       c->diag[i] = a->diag[i];
194817ab2063SBarry Smith     }
19493a40ed3dSBarry Smith   } else c->diag          = 0;
19506c7ebb05SLois Curfman McInnes   c->inode.limit        = a->inode.limit;
19516c7ebb05SLois Curfman McInnes   c->inode.max_limit    = a->inode.max_limit;
1952754ec7b1SSatish Balay   if (a->inode.size){
1953daed632aSSatish Balay     c->inode.size       = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size);
1954754ec7b1SSatish Balay     c->inode.node_count = a->inode.node_count;
1955daed632aSSatish Balay     PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int));
1956754ec7b1SSatish Balay   } else {
1957754ec7b1SSatish Balay     c->inode.size       = 0;
1958754ec7b1SSatish Balay     c->inode.node_count = 0;
1959754ec7b1SSatish Balay   }
1960416022c9SBarry Smith   c->nz                 = a->nz;
1961416022c9SBarry Smith   c->maxnz              = a->maxnz;
1962416022c9SBarry Smith   c->solve_work         = 0;
196376dd722bSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1964754ec7b1SSatish Balay 
1965416022c9SBarry Smith   *B = C;
19663a40ed3dSBarry Smith   PetscFunctionReturn(0);
196717ab2063SBarry Smith }
196817ab2063SBarry Smith 
19695615d1e5SSatish Balay #undef __FUNC__
19705615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ"
197119bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A)
197217ab2063SBarry Smith {
1973416022c9SBarry Smith   Mat_SeqAIJ   *a;
1974416022c9SBarry Smith   Mat          B;
197517699dbbSLois Curfman McInnes   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1976bcd2baecSBarry Smith   MPI_Comm     comm;
197717ab2063SBarry Smith 
19783a40ed3dSBarry Smith   PetscFunctionBegin;
197919bcc07fSBarry Smith   PetscObjectGetComm((PetscObject) viewer,&comm);
198017699dbbSLois Curfman McInnes   MPI_Comm_size(comm,&size);
1981a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor");
198290ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
19830752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
1984a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file");
198517ab2063SBarry Smith   M = header[1]; N = header[2]; nz = header[3];
198617ab2063SBarry Smith 
1987d64ed03dSBarry Smith   if (nz < 0) {
1988a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ");
1989d64ed03dSBarry Smith   }
1990d64ed03dSBarry Smith 
199117ab2063SBarry Smith   /* read in row lengths */
19920452661fSBarry Smith   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
19930752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
199417ab2063SBarry Smith 
199517ab2063SBarry Smith   /* create our matrix */
1996416022c9SBarry Smith   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1997416022c9SBarry Smith   B = *A;
1998416022c9SBarry Smith   a = (Mat_SeqAIJ *) B->data;
1999416022c9SBarry Smith   shift = a->indexshift;
200017ab2063SBarry Smith 
200117ab2063SBarry Smith   /* read in column indices and adjust for Fortran indexing*/
20020752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr);
200317ab2063SBarry Smith   if (shift) {
200417ab2063SBarry Smith     for ( i=0; i<nz; i++ ) {
2005416022c9SBarry Smith       a->j[i] += 1;
200617ab2063SBarry Smith     }
200717ab2063SBarry Smith   }
200817ab2063SBarry Smith 
200917ab2063SBarry Smith   /* read in nonzero values */
20100752156aSBarry Smith   ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr);
201117ab2063SBarry Smith 
201217ab2063SBarry Smith   /* set matrix "i" values */
2013416022c9SBarry Smith   a->i[0] = -shift;
201417ab2063SBarry Smith   for ( i=1; i<= M; i++ ) {
2015416022c9SBarry Smith     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2016416022c9SBarry Smith     a->ilen[i-1] = rowlengths[i-1];
201717ab2063SBarry Smith   }
20180452661fSBarry Smith   PetscFree(rowlengths);
201917ab2063SBarry Smith 
20206d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20216d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
20223a40ed3dSBarry Smith   PetscFunctionReturn(0);
202317ab2063SBarry Smith }
202417ab2063SBarry Smith 
20255615d1e5SSatish Balay #undef __FUNC__
20265615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ"
20278f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg)
20287264ac53SSatish Balay {
20297264ac53SSatish Balay   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data;
20307264ac53SSatish Balay 
20313a40ed3dSBarry Smith   PetscFunctionBegin;
2032a8c6a408SBarry Smith   if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
20337264ac53SSatish Balay 
20347264ac53SSatish Balay   /* If the  matrix dimensions are not equal, or no of nonzeros or shift */
20357264ac53SSatish Balay   if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)||
2036bcd2baecSBarry Smith       (a->indexshift != b->indexshift)) {
20373a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2038bcd2baecSBarry Smith   }
20397264ac53SSatish Balay 
20407264ac53SSatish Balay   /* if the a->i are the same */
20418108c231SLois Curfman McInnes   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
20423a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
20437264ac53SSatish Balay   }
20447264ac53SSatish Balay 
20457264ac53SSatish Balay   /* if a->j are the same */
2046bcd2baecSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
20473a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
2048bcd2baecSBarry Smith   }
2049bcd2baecSBarry Smith 
2050bcd2baecSBarry Smith   /* if a->a are the same */
205119bcc07fSBarry Smith   if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) {
20523a40ed3dSBarry Smith     *flg = PETSC_FALSE; PetscFunctionReturn(0);
20537264ac53SSatish Balay   }
205477c4ece6SBarry Smith   *flg = PETSC_TRUE;
20553a40ed3dSBarry Smith   PetscFunctionReturn(0);
20567264ac53SSatish Balay 
20577264ac53SSatish Balay }
2058