xref: /petsc/src/mat/impls/baij/seq/baij.c (revision f1e2ffcdf65679ed52caf95ee2fbac3db31d26d8)
1*f1e2ffcdSBarry Smith /*$Id: baij.c,v 1.192 1999/11/20 18:36:47 bsmith Exp bsmith $*/
22593348eSBarry Smith 
32593348eSBarry Smith /*
4b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
52593348eSBarry Smith   matrix storage format.
62593348eSBarry Smith */
73369ce9aSBarry Smith #include "sys.h"
870f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
91a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
101a6a6d98SBarry Smith #include "src/inline/spops.h"
113b547af2SSatish Balay 
122d61bbb3SSatish Balay #define CHUNKSIZE  10
13de6a44a3SBarry Smith 
14be5855fcSBarry Smith /*
15be5855fcSBarry Smith      Checks for missing diagonals
16be5855fcSBarry Smith */
17be5855fcSBarry Smith #undef __FUNC__
18c4992f7dSBarry Smith #define __FUNC__ "MatMissingDiagonal_SeqBAIJ"
19c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A)
20be5855fcSBarry Smith {
21be5855fcSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
22883fce79SBarry Smith   int         *diag, *jj = a->j,i,ierr;
23be5855fcSBarry Smith 
24be5855fcSBarry Smith   PetscFunctionBegin;
25c4992f7dSBarry Smith   ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
26883fce79SBarry Smith   diag = a->diag;
270e8e8aceSBarry Smith   for ( i=0; i<a->mbs; i++ ) {
28be5855fcSBarry Smith     if (jj[diag[i]] != i) {
29be5855fcSBarry Smith       SETERRQ1(1,1,"Matrix is missing diagonal number %d",i);
30be5855fcSBarry Smith     }
31be5855fcSBarry Smith   }
32be5855fcSBarry Smith   PetscFunctionReturn(0);
33be5855fcSBarry Smith }
34be5855fcSBarry Smith 
355615d1e5SSatish Balay #undef __FUNC__
36c4992f7dSBarry Smith #define __FUNC__ "MatMarkDiagonal_SeqBAIJ"
37c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A)
38de6a44a3SBarry Smith {
39de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
407fc0212eSBarry Smith   int         i,j, *diag, m = a->mbs;
41de6a44a3SBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
43883fce79SBarry Smith   if (a->diag) PetscFunctionReturn(0);
44883fce79SBarry Smith 
45de6a44a3SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag);
46de6a44a3SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
477fc0212eSBarry Smith   for ( i=0; i<m; i++ ) {
48ecc615b2SBarry Smith     diag[i] = a->i[i+1];
49de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
50de6a44a3SBarry Smith       if (a->j[j] == i) {
51de6a44a3SBarry Smith         diag[i] = j;
52de6a44a3SBarry Smith         break;
53de6a44a3SBarry Smith       }
54de6a44a3SBarry Smith     }
55de6a44a3SBarry Smith   }
56de6a44a3SBarry Smith   a->diag = diag;
573a40ed3dSBarry Smith   PetscFunctionReturn(0);
58de6a44a3SBarry Smith }
592593348eSBarry Smith 
602593348eSBarry Smith 
613b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
623b2fbd54SBarry Smith 
635615d1e5SSatish Balay #undef __FUNC__
645615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ"
653b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
663b2fbd54SBarry Smith                             PetscTruth *done)
673b2fbd54SBarry Smith {
683b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
693b2fbd54SBarry Smith   int         ierr, n = a->mbs,i;
703b2fbd54SBarry Smith 
713a40ed3dSBarry Smith   PetscFunctionBegin;
723b2fbd54SBarry Smith   *nn = n;
733a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
743b2fbd54SBarry Smith   if (symmetric) {
753b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
763b2fbd54SBarry Smith   } else if (oshift == 1) {
773b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
783b2fbd54SBarry Smith     int nz = a->i[n] + 1;
793b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]++;
803b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]++;
813b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
823b2fbd54SBarry Smith   } else {
833b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
843b2fbd54SBarry Smith   }
853b2fbd54SBarry Smith 
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
873b2fbd54SBarry Smith }
883b2fbd54SBarry Smith 
895615d1e5SSatish Balay #undef __FUNC__
90d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ"
913b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
923b2fbd54SBarry Smith                                 PetscTruth *done)
933b2fbd54SBarry Smith {
943b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
95606d414cSSatish Balay   int         i,n = a->mbs,ierr;
963b2fbd54SBarry Smith 
973a40ed3dSBarry Smith   PetscFunctionBegin;
983a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
993b2fbd54SBarry Smith   if (symmetric) {
100606d414cSSatish Balay     ierr = PetscFree(*ia);CHKERRQ(ierr);
101606d414cSSatish Balay     ierr = PetscFree(*ja);CHKERRQ(ierr);
102af5da2bfSBarry Smith   } else if (oshift == 1) {
1033b2fbd54SBarry Smith     int nz = a->i[n];
1043b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]--;
1053b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]--;
1063b2fbd54SBarry Smith   }
1073a40ed3dSBarry Smith   PetscFunctionReturn(0);
1083b2fbd54SBarry Smith }
1093b2fbd54SBarry Smith 
1102d61bbb3SSatish Balay #undef __FUNC__
1112d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ"
1122d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs)
1132d61bbb3SSatish Balay {
1142d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data;
1152d61bbb3SSatish Balay 
1162d61bbb3SSatish Balay   PetscFunctionBegin;
1172d61bbb3SSatish Balay   *bs = baij->bs;
1182d61bbb3SSatish Balay   PetscFunctionReturn(0);
1192d61bbb3SSatish Balay }
1202d61bbb3SSatish Balay 
1212d61bbb3SSatish Balay 
1222d61bbb3SSatish Balay #undef __FUNC__
1232d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ"
124e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
1252d61bbb3SSatish Balay {
1262d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
127e51c0b9cSSatish Balay   int         ierr;
1282d61bbb3SSatish Balay 
129433994e6SBarry Smith   PetscFunctionBegin;
13094d884c6SBarry Smith   if (A->mapping) {
13194d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr);
13294d884c6SBarry Smith   }
13394d884c6SBarry Smith   if (A->bmapping) {
13494d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr);
13594d884c6SBarry Smith   }
13661b13de0SBarry Smith   if (A->rmap) {
13761b13de0SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
13861b13de0SBarry Smith   }
13961b13de0SBarry Smith   if (A->cmap) {
14061b13de0SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
14161b13de0SBarry Smith   }
142aa482453SBarry Smith #if defined(PETSC_USE_LOG)
143e1311b90SBarry Smith   PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
1442d61bbb3SSatish Balay #endif
145606d414cSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
146606d414cSSatish Balay   if (!a->singlemalloc) {
147606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
148606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
149606d414cSSatish Balay   }
150c38d4ed2SBarry Smith   if (a->row) {
151c38d4ed2SBarry Smith     ierr = ISDestroy(a->row);CHKERRQ(ierr);
152c38d4ed2SBarry Smith   }
153c38d4ed2SBarry Smith   if (a->col) {
154c38d4ed2SBarry Smith     ierr = ISDestroy(a->col);CHKERRQ(ierr);
155c38d4ed2SBarry Smith   }
156606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
157606d414cSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
158606d414cSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
159606d414cSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
160606d414cSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
161e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
162606d414cSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
163606d414cSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
1642d61bbb3SSatish Balay   PLogObjectDestroy(A);
1652d61bbb3SSatish Balay   PetscHeaderDestroy(A);
1662d61bbb3SSatish Balay   PetscFunctionReturn(0);
1672d61bbb3SSatish Balay }
1682d61bbb3SSatish Balay 
1692d61bbb3SSatish Balay #undef __FUNC__
1702d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ"
1712d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1722d61bbb3SSatish Balay {
1732d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1742d61bbb3SSatish Balay 
1752d61bbb3SSatish Balay   PetscFunctionBegin;
176c4992f7dSBarry Smith   if      (op == MAT_ROW_ORIENTED)                 a->roworiented    = PETSC_TRUE;
177c4992f7dSBarry Smith   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented    = PETSC_FALSE;
178c4992f7dSBarry Smith   else if (op == MAT_COLUMNS_SORTED)               a->sorted         = PETSC_TRUE;
179c4992f7dSBarry Smith   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted         = PETSC_FALSE;
180c4992f7dSBarry Smith   else if (op == MAT_KEEP_ZEROED_ROWS)             a->keepzeroedrows = PETSC_TRUE;
1812d61bbb3SSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew          = 1;
1824787f768SSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew          = -1;
1834787f768SSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew          = -2;
1842d61bbb3SSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew          = 0;
1852d61bbb3SSatish Balay   else if (op == MAT_ROWS_SORTED ||
1862d61bbb3SSatish Balay            op == MAT_ROWS_UNSORTED ||
1872d61bbb3SSatish Balay            op == MAT_SYMMETRIC ||
1882d61bbb3SSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
1892d61bbb3SSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
190b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
191b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE) {
192981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
1932d61bbb3SSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
1942d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
1952d61bbb3SSatish Balay   } else {
1962d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
1972d61bbb3SSatish Balay   }
1982d61bbb3SSatish Balay   PetscFunctionReturn(0);
1992d61bbb3SSatish Balay }
2002d61bbb3SSatish Balay 
2012d61bbb3SSatish Balay 
2022d61bbb3SSatish Balay #undef __FUNC__
2032d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ"
2042d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
2052d61bbb3SSatish Balay {
2062d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2072d61bbb3SSatish Balay 
2082d61bbb3SSatish Balay   PetscFunctionBegin;
2092d61bbb3SSatish Balay   if (m) *m = a->m;
2102d61bbb3SSatish Balay   if (n) *n = a->n;
2112d61bbb3SSatish Balay   PetscFunctionReturn(0);
2122d61bbb3SSatish Balay }
2132d61bbb3SSatish Balay 
2142d61bbb3SSatish Balay #undef __FUNC__
2152d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ"
2162d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
2172d61bbb3SSatish Balay {
2182d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2192d61bbb3SSatish Balay 
2202d61bbb3SSatish Balay   PetscFunctionBegin;
2212d61bbb3SSatish Balay   *m = 0; *n = a->m;
2222d61bbb3SSatish Balay   PetscFunctionReturn(0);
2232d61bbb3SSatish Balay }
2242d61bbb3SSatish Balay 
2252d61bbb3SSatish Balay #undef __FUNC__
2262d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ"
2272d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2282d61bbb3SSatish Balay {
2292d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *) A->data;
2302d61bbb3SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
2313f1db9ecSBarry Smith   MatScalar    *aa,*aa_i;
2323f1db9ecSBarry Smith   Scalar       *v_i;
2332d61bbb3SSatish Balay 
2342d61bbb3SSatish Balay   PetscFunctionBegin;
2352d61bbb3SSatish Balay   bs  = a->bs;
2362d61bbb3SSatish Balay   ai  = a->i;
2372d61bbb3SSatish Balay   aj  = a->j;
2382d61bbb3SSatish Balay   aa  = a->a;
2392d61bbb3SSatish Balay   bs2 = a->bs2;
2402d61bbb3SSatish Balay 
2412d61bbb3SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
2422d61bbb3SSatish Balay 
2432d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2442d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2452d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2462d61bbb3SSatish Balay   *nz = bs*M;
2472d61bbb3SSatish Balay 
2482d61bbb3SSatish Balay   if (v) {
2492d61bbb3SSatish Balay     *v = 0;
2502d61bbb3SSatish Balay     if (*nz) {
2512d61bbb3SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) );CHKPTRQ(*v);
2522d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2532d61bbb3SSatish Balay         v_i  = *v + i*bs;
2542d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2552d61bbb3SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
2562d61bbb3SSatish Balay       }
2572d61bbb3SSatish Balay     }
2582d61bbb3SSatish Balay   }
2592d61bbb3SSatish Balay 
2602d61bbb3SSatish Balay   if (idx) {
2612d61bbb3SSatish Balay     *idx = 0;
2622d61bbb3SSatish Balay     if (*nz) {
2632d61bbb3SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) );CHKPTRQ(*idx);
2642d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2652d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2662d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2672d61bbb3SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
2682d61bbb3SSatish Balay       }
2692d61bbb3SSatish Balay     }
2702d61bbb3SSatish Balay   }
2712d61bbb3SSatish Balay   PetscFunctionReturn(0);
2722d61bbb3SSatish Balay }
2732d61bbb3SSatish Balay 
2742d61bbb3SSatish Balay #undef __FUNC__
2752d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ"
2762d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2772d61bbb3SSatish Balay {
278606d414cSSatish Balay   int ierr;
279606d414cSSatish Balay 
2802d61bbb3SSatish Balay   PetscFunctionBegin;
281606d414cSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
282606d414cSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
2832d61bbb3SSatish Balay   PetscFunctionReturn(0);
2842d61bbb3SSatish Balay }
2852d61bbb3SSatish Balay 
2862d61bbb3SSatish Balay #undef __FUNC__
2872d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ"
2882d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
2892d61bbb3SSatish Balay {
2902d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
2912d61bbb3SSatish Balay   Mat         C;
2922d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
2932d61bbb3SSatish Balay   int         *rows,*cols,bs2=a->bs2;
2943f1db9ecSBarry Smith   MatScalar   *array = a->a;
2952d61bbb3SSatish Balay 
2962d61bbb3SSatish Balay   PetscFunctionBegin;
2972d61bbb3SSatish Balay   if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
2982d61bbb3SSatish Balay   col  = (int *) PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col);
299549d3d68SSatish Balay   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
3002d61bbb3SSatish Balay 
3012d61bbb3SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
3022d61bbb3SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
303606d414cSSatish Balay   ierr = PetscFree(col);CHKERRQ(ierr);
3042d61bbb3SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows);
3052d61bbb3SSatish Balay   cols = rows + bs;
3062d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3072d61bbb3SSatish Balay     cols[0] = i*bs;
3082d61bbb3SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
3092d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
3102d61bbb3SSatish Balay     for ( j=0; j<len; j++ ) {
3112d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
3122d61bbb3SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
3132d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
3142d61bbb3SSatish Balay       array += bs2;
3152d61bbb3SSatish Balay     }
3162d61bbb3SSatish Balay   }
317606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
3182d61bbb3SSatish Balay 
3192d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3202d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3212d61bbb3SSatish Balay 
322c4992f7dSBarry Smith   if (B) {
3232d61bbb3SSatish Balay     *B = C;
3242d61bbb3SSatish Balay   } else {
325f830108cSBarry Smith     PetscOps *Abops;
326cc2dc46cSBarry Smith     MatOps   Aops;
327f830108cSBarry Smith 
3282d61bbb3SSatish Balay     /* This isn't really an in-place transpose */
329606d414cSSatish Balay     ierr = PetscFree(a->a);CHKERRQ(ierr);
330606d414cSSatish Balay     if (!a->singlemalloc) {
331606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
332606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
333606d414cSSatish Balay     }
334606d414cSSatish Balay     if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
335606d414cSSatish Balay     if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
336606d414cSSatish Balay     if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
337606d414cSSatish Balay     if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
338606d414cSSatish Balay     ierr = PetscFree(a);CHKERRQ(ierr);
339f830108cSBarry Smith 
3407413bad6SBarry Smith 
3417413bad6SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
3427413bad6SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
3437413bad6SBarry Smith 
344f830108cSBarry Smith     /*
345f830108cSBarry Smith        This is horrible, horrible code. We need to keep the
346f830108cSBarry Smith       A pointers for the bops and ops but copy everything
347f830108cSBarry Smith       else from C.
348f830108cSBarry Smith     */
349f830108cSBarry Smith     Abops    = A->bops;
350f830108cSBarry Smith     Aops     = A->ops;
351549d3d68SSatish Balay     ierr     = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
352f830108cSBarry Smith     A->bops  = Abops;
353f830108cSBarry Smith     A->ops   = Aops;
35427a8da17SBarry Smith     A->qlist = 0;
355f830108cSBarry Smith 
3562d61bbb3SSatish Balay     PetscHeaderDestroy(C);
3572d61bbb3SSatish Balay   }
3582d61bbb3SSatish Balay   PetscFunctionReturn(0);
3592d61bbb3SSatish Balay }
3602d61bbb3SSatish Balay 
3615615d1e5SSatish Balay #undef __FUNC__
362d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary"
363b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
3642593348eSBarry Smith {
365b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3669df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
367b6490206SBarry Smith   Scalar      *aa;
368ce6f0cecSBarry Smith   FILE        *file;
3692593348eSBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
37190ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
3722593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
3732593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3743b2fbd54SBarry Smith 
3752593348eSBarry Smith   col_lens[1] = a->m;
3762593348eSBarry Smith   col_lens[2] = a->n;
3777e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3782593348eSBarry Smith 
3792593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
380b6490206SBarry Smith   count = 0;
381b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
382b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
383b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
384b6490206SBarry Smith     }
3852593348eSBarry Smith   }
3860752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr);
387606d414cSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
3882593348eSBarry Smith 
3892593348eSBarry Smith   /* store column indices (zero start index) */
39066963ce1SSatish Balay   jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) );CHKPTRQ(jj);
391b6490206SBarry Smith   count = 0;
392b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
393b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
394b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
395b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
396b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3972593348eSBarry Smith         }
3982593348eSBarry Smith       }
399b6490206SBarry Smith     }
400b6490206SBarry Smith   }
4010752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
402606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
4032593348eSBarry Smith 
4042593348eSBarry Smith   /* store nonzero values */
40566963ce1SSatish Balay   aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa);
406b6490206SBarry Smith   count = 0;
407b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
408b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
409b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
410b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
4117e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
412b6490206SBarry Smith         }
413b6490206SBarry Smith       }
414b6490206SBarry Smith     }
415b6490206SBarry Smith   }
4160752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
417606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
418ce6f0cecSBarry Smith 
419ce6f0cecSBarry Smith   ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
420ce6f0cecSBarry Smith   if (file) {
421ce6f0cecSBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
422ce6f0cecSBarry Smith   }
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
4242593348eSBarry Smith }
4252593348eSBarry Smith 
4265615d1e5SSatish Balay #undef __FUNC__
427d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII"
428b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
4292593348eSBarry Smith {
430b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4319df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
4322593348eSBarry Smith   char        *outputname;
4332593348eSBarry Smith 
4343a40ed3dSBarry Smith   PetscFunctionBegin;
43577ed5343SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
436888f2ed8SSatish Balay   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
437639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
438d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
4390ef38995SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
4406831982aSBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
4410ef38995SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
4426831982aSBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
44344cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
44444cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
4456831982aSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
44644cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
44744cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
448aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4490ef38995SBarry Smith             if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
4506831982aSBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
4516831982aSBarry Smith                       PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4520ef38995SBarry Smith             } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
4536831982aSBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
4546831982aSBarry Smith                       PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4550ef38995SBarry Smith             } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
4566831982aSBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4570ef38995SBarry Smith             }
45844cd7ae7SLois Curfman McInnes #else
4590ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
4606831982aSBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
4610ef38995SBarry Smith             }
46244cd7ae7SLois Curfman McInnes #endif
46344cd7ae7SLois Curfman McInnes           }
46444cd7ae7SLois Curfman McInnes         }
4656831982aSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
46644cd7ae7SLois Curfman McInnes       }
46744cd7ae7SLois Curfman McInnes     }
4686831982aSBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
4690ef38995SBarry Smith   } else {
4706831982aSBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
471b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
472b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
4736831982aSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
474b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
475b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
476aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
477e20fef11SSatish Balay             if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) {
4786831982aSBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
4796831982aSBarry Smith                 PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4800ef38995SBarry Smith             } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) {
4816831982aSBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
4826831982aSBarry Smith                 PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4830ef38995SBarry Smith             } else {
4846831982aSBarry Smith               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
48588685aaeSLois Curfman McInnes             }
48688685aaeSLois Curfman McInnes #else
4876831982aSBarry Smith             ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
48888685aaeSLois Curfman McInnes #endif
4892593348eSBarry Smith           }
4902593348eSBarry Smith         }
4916831982aSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
4922593348eSBarry Smith       }
4932593348eSBarry Smith     }
4946831982aSBarry Smith     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
495b6490206SBarry Smith   }
4966831982aSBarry Smith   ierr = ViewerFlush(viewer);CHKERRQ(ierr);
4973a40ed3dSBarry Smith   PetscFunctionReturn(0);
4982593348eSBarry Smith }
4992593348eSBarry Smith 
5005615d1e5SSatish Balay #undef __FUNC__
50177ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom"
50277ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa)
5033270192aSSatish Balay {
50477ed5343SBarry Smith   Mat          A = (Mat) Aa;
5053270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
5067dce120fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
507fae8c767SSatish Balay   double       xl,yl,xr,yr,x_l,x_r,y_l,y_r;
5083f1db9ecSBarry Smith   MatScalar    *aa;
50977ed5343SBarry Smith   MPI_Comm     comm;
51077ed5343SBarry Smith   Viewer       viewer;
5113270192aSSatish Balay 
5123a40ed3dSBarry Smith   PetscFunctionBegin;
51377ed5343SBarry Smith   /*
51477ed5343SBarry Smith       This is nasty. If this is called from an originally parallel matrix
51577ed5343SBarry Smith    then all processes call this, but only the first has the matrix so the
51677ed5343SBarry Smith    rest should return immediately.
51777ed5343SBarry Smith   */
51877ed5343SBarry Smith   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
519d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
52077ed5343SBarry Smith   if (rank) PetscFunctionReturn(0);
5213270192aSSatish Balay 
52277ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr);
52377ed5343SBarry Smith 
52477ed5343SBarry Smith   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
52577ed5343SBarry Smith 
5263270192aSSatish Balay   /* loop over matrix elements drawing boxes */
5273270192aSSatish Balay   color = DRAW_BLUE;
5283270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5293270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5303270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5313270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5323270192aSSatish Balay       aa = a->a + j*bs2;
5333270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5343270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5353270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
536433994e6SBarry Smith           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5373270192aSSatish Balay         }
5383270192aSSatish Balay       }
5393270192aSSatish Balay     }
5403270192aSSatish Balay   }
5413270192aSSatish Balay   color = DRAW_CYAN;
5423270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5433270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5443270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5453270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5463270192aSSatish Balay       aa = a->a + j*bs2;
5473270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5483270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5493270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
550433994e6SBarry Smith           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5513270192aSSatish Balay         }
5523270192aSSatish Balay       }
5533270192aSSatish Balay     }
5543270192aSSatish Balay   }
5553270192aSSatish Balay 
5563270192aSSatish Balay   color = DRAW_RED;
5573270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5583270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5593270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5603270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5613270192aSSatish Balay       aa = a->a + j*bs2;
5623270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5633270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5643270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
565433994e6SBarry Smith           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5663270192aSSatish Balay         }
5673270192aSSatish Balay       }
5683270192aSSatish Balay     }
5693270192aSSatish Balay   }
57077ed5343SBarry Smith   PetscFunctionReturn(0);
57177ed5343SBarry Smith }
5723270192aSSatish Balay 
57377ed5343SBarry Smith #undef __FUNC__
57477ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw"
57577ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
57677ed5343SBarry Smith {
57777ed5343SBarry Smith   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
5787dce120fSSatish Balay   int          ierr;
5797dce120fSSatish Balay   double       xl,yl,xr,yr,w,h;
58077ed5343SBarry Smith   Draw         draw;
58177ed5343SBarry Smith   PetscTruth   isnull;
5823270192aSSatish Balay 
58377ed5343SBarry Smith   PetscFunctionBegin;
58477ed5343SBarry Smith 
58577ed5343SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
58677ed5343SBarry Smith   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
58777ed5343SBarry Smith 
58877ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
58977ed5343SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
59077ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
5913270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
59277ed5343SBarry Smith   ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
59377ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5943a40ed3dSBarry Smith   PetscFunctionReturn(0);
5953270192aSSatish Balay }
5963270192aSSatish Balay 
5975615d1e5SSatish Balay #undef __FUNC__
598d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ"
599e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer)
6002593348eSBarry Smith {
60119bcc07fSBarry Smith   int        ierr;
6026831982aSBarry Smith   PetscTruth issocket,isascii,isbinary,isdraw;
6032593348eSBarry Smith 
6043a40ed3dSBarry Smith   PetscFunctionBegin;
6056831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
6066831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
6076831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
6086831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
6090f5bd95cSBarry Smith   if (issocket) {
6107b2a1423SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported");
6110f5bd95cSBarry Smith   } else if (isascii){
6123a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
6130f5bd95cSBarry Smith   } else if (isbinary) {
6143a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
6150f5bd95cSBarry Smith   } else if (isdraw) {
6163a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
6175cd90555SBarry Smith   } else {
6180f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
6192593348eSBarry Smith   }
6203a40ed3dSBarry Smith   PetscFunctionReturn(0);
6212593348eSBarry Smith }
622b6490206SBarry Smith 
623cd0e1443SSatish Balay 
6245615d1e5SSatish Balay #undef __FUNC__
6252d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ"
6262d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
627cd0e1443SSatish Balay {
628cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6292d61bbb3SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
6302d61bbb3SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
6312d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6323f1db9ecSBarry Smith   MatScalar  *ap, *aa = a->a, zero = 0.0;
633cd0e1443SSatish Balay 
6343a40ed3dSBarry Smith   PetscFunctionBegin;
6352d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
636cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
637a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
638a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
6392d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6402c3acbe9SBarry Smith     nrow = ailen[brow];
6412d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
642a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
643a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6442d61bbb3SSatish Balay       col  = in[l] ;
6452d61bbb3SSatish Balay       bcol = col/bs;
6462d61bbb3SSatish Balay       cidx = col%bs;
6472d61bbb3SSatish Balay       ridx = row%bs;
6482d61bbb3SSatish Balay       high = nrow;
6492d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6502d61bbb3SSatish Balay       while (high-low > 5) {
651cd0e1443SSatish Balay         t = (low+high)/2;
652cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
653cd0e1443SSatish Balay         else             low  = t;
654cd0e1443SSatish Balay       }
655cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
656cd0e1443SSatish Balay         if (rp[i] > bcol) break;
657cd0e1443SSatish Balay         if (rp[i] == bcol) {
6582d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6592d61bbb3SSatish Balay           goto finished;
660cd0e1443SSatish Balay         }
661cd0e1443SSatish Balay       }
6622d61bbb3SSatish Balay       *v++ = zero;
6632d61bbb3SSatish Balay       finished:;
664cd0e1443SSatish Balay     }
665cd0e1443SSatish Balay   }
6663a40ed3dSBarry Smith   PetscFunctionReturn(0);
667cd0e1443SSatish Balay }
668cd0e1443SSatish Balay 
6692d61bbb3SSatish Balay 
6705615d1e5SSatish Balay #undef __FUNC__
67105a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ"
67292c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
67392c4ed94SBarry Smith {
67492c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6758a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
676dd9472c6SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented;
677549d3d68SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
6783f1db9ecSBarry Smith   Scalar      *value = v;
6793f1db9ecSBarry Smith   MatScalar   *ap,*aa=a->a,*bap;
68092c4ed94SBarry Smith 
6813a40ed3dSBarry Smith   PetscFunctionBegin;
6820e324ae4SSatish Balay   if (roworiented) {
6830e324ae4SSatish Balay     stepval = (n-1)*bs;
6840e324ae4SSatish Balay   } else {
6850e324ae4SSatish Balay     stepval = (m-1)*bs;
6860e324ae4SSatish Balay   }
68792c4ed94SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
68892c4ed94SBarry Smith     row  = im[k];
6895ef9f2a5SBarry Smith     if (row < 0) continue;
690aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
691a8c6a408SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
69292c4ed94SBarry Smith #endif
69392c4ed94SBarry Smith     rp   = aj + ai[row];
69492c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
69592c4ed94SBarry Smith     rmax = imax[row];
69692c4ed94SBarry Smith     nrow = ailen[row];
69792c4ed94SBarry Smith     low  = 0;
69892c4ed94SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
6995ef9f2a5SBarry Smith       if (in[l] < 0) continue;
700aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
701a8c6a408SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
70292c4ed94SBarry Smith #endif
70392c4ed94SBarry Smith       col = in[l];
70492c4ed94SBarry Smith       if (roworiented) {
7050e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
7060e324ae4SSatish Balay       } else {
7070e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
70892c4ed94SBarry Smith       }
70992c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
71092c4ed94SBarry Smith       while (high-low > 7) {
71192c4ed94SBarry Smith         t = (low+high)/2;
71292c4ed94SBarry Smith         if (rp[t] > col) high = t;
71392c4ed94SBarry Smith         else             low  = t;
71492c4ed94SBarry Smith       }
71592c4ed94SBarry Smith       for ( i=low; i<high; i++ ) {
71692c4ed94SBarry Smith         if (rp[i] > col) break;
71792c4ed94SBarry Smith         if (rp[i] == col) {
7188a84c255SSatish Balay           bap  = ap +  bs2*i;
7190e324ae4SSatish Balay           if (roworiented) {
7208a84c255SSatish Balay             if (is == ADD_VALUES) {
721dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
722dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7238a84c255SSatish Balay                   bap[jj] += *value++;
724dd9472c6SBarry Smith                 }
725dd9472c6SBarry Smith               }
7260e324ae4SSatish Balay             } else {
727dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
728dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7290e324ae4SSatish Balay                   bap[jj] = *value++;
7308a84c255SSatish Balay                 }
731dd9472c6SBarry Smith               }
732dd9472c6SBarry Smith             }
7330e324ae4SSatish Balay           } else {
7340e324ae4SSatish Balay             if (is == ADD_VALUES) {
735dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
736dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7370e324ae4SSatish Balay                   *bap++ += *value++;
738dd9472c6SBarry Smith                 }
739dd9472c6SBarry Smith               }
7400e324ae4SSatish Balay             } else {
741dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
742dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7430e324ae4SSatish Balay                   *bap++  = *value++;
7440e324ae4SSatish Balay                 }
7458a84c255SSatish Balay               }
746dd9472c6SBarry Smith             }
747dd9472c6SBarry Smith           }
748f1241b54SBarry Smith           goto noinsert2;
74992c4ed94SBarry Smith         }
75092c4ed94SBarry Smith       }
75189280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
752a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
75392c4ed94SBarry Smith       if (nrow >= rmax) {
75492c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
75592c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
7563f1db9ecSBarry Smith         MatScalar *new_a;
75792c4ed94SBarry Smith 
758a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
75989280ab3SLois Curfman McInnes 
76092c4ed94SBarry Smith         /* malloc new storage space */
7613f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
7623f1db9ecSBarry Smith         new_a   = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a);
76392c4ed94SBarry Smith         new_j   = (int *) (new_a + bs2*new_nz);
76492c4ed94SBarry Smith         new_i   = new_j + new_nz;
76592c4ed94SBarry Smith 
76692c4ed94SBarry Smith         /* copy over old data into new slots */
76792c4ed94SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
76892c4ed94SBarry Smith         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
769549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
77092c4ed94SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
771549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
772549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
773549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
774549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
77592c4ed94SBarry Smith         /* free up old matrix storage */
776606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
777606d414cSSatish Balay         if (!a->singlemalloc) {
778606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
779606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
780606d414cSSatish Balay         }
78192c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
782c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
78392c4ed94SBarry Smith 
78492c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
78592c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
7863f1db9ecSBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
78792c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
78892c4ed94SBarry Smith         a->reallocs++;
78992c4ed94SBarry Smith         a->nz++;
79092c4ed94SBarry Smith       }
79192c4ed94SBarry Smith       N = nrow++ - 1;
79292c4ed94SBarry Smith       /* shift up all the later entries in this row */
79392c4ed94SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
79492c4ed94SBarry Smith         rp[ii+1] = rp[ii];
795549d3d68SSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
79692c4ed94SBarry Smith       }
797549d3d68SSatish Balay       if (N >= i) {
798549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
799549d3d68SSatish Balay       }
80092c4ed94SBarry Smith       rp[i] = col;
8018a84c255SSatish Balay       bap   = ap +  bs2*i;
8020e324ae4SSatish Balay       if (roworiented) {
803dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval) {
804dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs ) {
8050e324ae4SSatish Balay             bap[jj] = *value++;
806dd9472c6SBarry Smith           }
807dd9472c6SBarry Smith         }
8080e324ae4SSatish Balay       } else {
809dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval ) {
810dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++ ) {
8110e324ae4SSatish Balay             *bap++  = *value++;
8120e324ae4SSatish Balay           }
813dd9472c6SBarry Smith         }
814dd9472c6SBarry Smith       }
815f1241b54SBarry Smith       noinsert2:;
81692c4ed94SBarry Smith       low = i;
81792c4ed94SBarry Smith     }
81892c4ed94SBarry Smith     ailen[row] = nrow;
81992c4ed94SBarry Smith   }
8203a40ed3dSBarry Smith   PetscFunctionReturn(0);
82192c4ed94SBarry Smith }
82292c4ed94SBarry Smith 
823f2501298SSatish Balay 
8245615d1e5SSatish Balay #undef __FUNC__
8255615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ"
8268f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
827584200bdSSatish Balay {
828584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
829584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
830a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
831549d3d68SSatish Balay   int        mbs = a->mbs, bs2 = a->bs2,rmax = 0,ierr;
8323f1db9ecSBarry Smith   MatScalar  *aa = a->a, *ap;
833584200bdSSatish Balay 
8343a40ed3dSBarry Smith   PetscFunctionBegin;
8353a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
836584200bdSSatish Balay 
83743ee02c3SBarry Smith   if (m) rmax = ailen[0];
838584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
839584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
840584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
841d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
842584200bdSSatish Balay     if (fshift) {
843a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
844584200bdSSatish Balay       N = ailen[i];
845584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
846584200bdSSatish Balay         ip[j-fshift] = ip[j];
847549d3d68SSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
848584200bdSSatish Balay       }
849584200bdSSatish Balay     }
850584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
851584200bdSSatish Balay   }
852584200bdSSatish Balay   if (mbs) {
853584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
854584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
855584200bdSSatish Balay   }
856584200bdSSatish Balay   /* reset ilen and imax for each row */
857584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
858584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
859584200bdSSatish Balay   }
860a7c10996SSatish Balay   a->nz = ai[mbs];
861584200bdSSatish Balay 
862584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
863584200bdSSatish Balay   if (fshift && a->diag) {
864606d414cSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
865584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
866584200bdSSatish Balay     a->diag = 0;
867584200bdSSatish Balay   }
8683d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
869219d9a1aSLois Curfman McInnes            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
8703d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
871584200bdSSatish Balay            a->reallocs);
872d402145bSBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
873e2f3b5e9SSatish Balay   a->reallocs          = 0;
8744e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
8754e220ebcSLois Curfman McInnes 
8763a40ed3dSBarry Smith   PetscFunctionReturn(0);
877584200bdSSatish Balay }
878584200bdSSatish Balay 
8795a838052SSatish Balay 
880bea157c4SSatish Balay 
881bea157c4SSatish Balay /*
882bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
883bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
884bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
885bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
886bea157c4SSatish Balay */
8875615d1e5SSatish Balay #undef __FUNC__
888bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks"
889bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
890d9b7c43dSSatish Balay {
891bea157c4SSatish Balay   int        i,j,k,row;
892bea157c4SSatish Balay   PetscTruth flg;
8933a40ed3dSBarry Smith 
894433994e6SBarry Smith   PetscFunctionBegin;
895bea157c4SSatish Balay   for ( i=0,j=0; i<n; j++ ) {
896bea157c4SSatish Balay     row = idx[i];
897bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
898bea157c4SSatish Balay       sizes[j] = 1;
899bea157c4SSatish Balay       i++;
900e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
901bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
902bea157c4SSatish Balay       i++;
903bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
904bea157c4SSatish Balay       flg = PETSC_TRUE;
905bea157c4SSatish Balay       for ( k=1; k<bs; k++ ) {
906bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
907bea157c4SSatish Balay           flg = PETSC_FALSE;
908bea157c4SSatish Balay           break;
909d9b7c43dSSatish Balay         }
910bea157c4SSatish Balay       }
911bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
912bea157c4SSatish Balay         sizes[j] = bs;
913bea157c4SSatish Balay         i+= bs;
914bea157c4SSatish Balay       } else {
915bea157c4SSatish Balay         sizes[j] = 1;
916bea157c4SSatish Balay         i++;
917bea157c4SSatish Balay       }
918bea157c4SSatish Balay     }
919bea157c4SSatish Balay   }
920bea157c4SSatish Balay   *bs_max = j;
9213a40ed3dSBarry Smith   PetscFunctionReturn(0);
922d9b7c43dSSatish Balay }
923d9b7c43dSSatish Balay 
9245615d1e5SSatish Balay #undef __FUNC__
9255615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ"
9268f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
927d9b7c43dSSatish Balay {
928d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
929b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
930bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
9313f1db9ecSBarry Smith   Scalar      zero = 0.0;
9323f1db9ecSBarry Smith   MatScalar   *aa;
933d9b7c43dSSatish Balay 
9343a40ed3dSBarry Smith   PetscFunctionBegin;
935d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
936d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
937d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
938d9b7c43dSSatish Balay 
939bea157c4SSatish Balay   /* allocate memory for rows,sizes */
940bea157c4SSatish Balay   rows  = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows);
941bea157c4SSatish Balay   sizes = rows + is_n;
942bea157c4SSatish Balay 
943bea157c4SSatish Balay   /* initialize copy IS valurs to rows, and sort them */
944bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
945bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
946dffd3267SBarry Smith   if (baij->keepzeroedrows) {
947dffd3267SBarry Smith     for ( i=0; i<is_n; i++ ) { sizes[i] = 1; }
948dffd3267SBarry Smith     bs_max = is_n;
949dffd3267SBarry Smith   } else {
950bea157c4SSatish Balay     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
951dffd3267SBarry Smith   }
952b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
953bea157c4SSatish Balay 
954bea157c4SSatish Balay   for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) {
955bea157c4SSatish Balay     row   = rows[j];
956b6815fffSBarry Smith     if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row);
957bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
958bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
959dffd3267SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
960bea157c4SSatish Balay       if (diag) {
961bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
962bea157c4SSatish Balay           baij->ilen[row/bs] = 1;
963bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
964bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
965a07cd24cSSatish Balay         }
966bea157c4SSatish Balay         /* Now insert all the diagoanl values for this bs */
967bea157c4SSatish Balay         for ( k=0; k<bs; k++ ) {
968bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
969bea157c4SSatish Balay         }
970bea157c4SSatish Balay       } else { /* (!diag) */
971bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
972bea157c4SSatish Balay       } /* end (!diag) */
973bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
974aa482453SBarry Smith #if defined (PETSC_USE_DEBUG)
975bea157c4SSatish Balay       if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1");
976bea157c4SSatish Balay #endif
977bea157c4SSatish Balay       for ( k=0; k<count; k++ ) {
978d9b7c43dSSatish Balay         aa[0] = zero;
979d9b7c43dSSatish Balay         aa+=bs;
980d9b7c43dSSatish Balay       }
981d9b7c43dSSatish Balay       if (diag) {
982f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
983d9b7c43dSSatish Balay       }
984d9b7c43dSSatish Balay     }
985bea157c4SSatish Balay   }
986bea157c4SSatish Balay 
987606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
9889a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9893a40ed3dSBarry Smith   PetscFunctionReturn(0);
990d9b7c43dSSatish Balay }
9911c351548SSatish Balay 
9925615d1e5SSatish Balay #undef __FUNC__
9932d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ"
9942d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
9952d61bbb3SSatish Balay {
9962d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9972d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
9982d61bbb3SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
9992d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1000549d3d68SSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
10013f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
10022d61bbb3SSatish Balay 
10032d61bbb3SSatish Balay   PetscFunctionBegin;
10042d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
10052d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
10065ef9f2a5SBarry Smith     if (row < 0) continue;
1007aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
100832e87ba3SBarry Smith     if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m);
10092d61bbb3SSatish Balay #endif
10102d61bbb3SSatish Balay     rp   = aj + ai[brow];
10112d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
10122d61bbb3SSatish Balay     rmax = imax[brow];
10132d61bbb3SSatish Balay     nrow = ailen[brow];
10142d61bbb3SSatish Balay     low  = 0;
10152d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
10165ef9f2a5SBarry Smith       if (in[l] < 0) continue;
1017aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
101832e87ba3SBarry Smith       if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n);
10192d61bbb3SSatish Balay #endif
10202d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
10212d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
10222d61bbb3SSatish Balay       if (roworiented) {
10235ef9f2a5SBarry Smith         value = v[l + k*n];
10242d61bbb3SSatish Balay       } else {
10252d61bbb3SSatish Balay         value = v[k + l*m];
10262d61bbb3SSatish Balay       }
10272d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
10282d61bbb3SSatish Balay       while (high-low > 7) {
10292d61bbb3SSatish Balay         t = (low+high)/2;
10302d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
10312d61bbb3SSatish Balay         else              low  = t;
10322d61bbb3SSatish Balay       }
10332d61bbb3SSatish Balay       for ( i=low; i<high; i++ ) {
10342d61bbb3SSatish Balay         if (rp[i] > bcol) break;
10352d61bbb3SSatish Balay         if (rp[i] == bcol) {
10362d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
10372d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
10382d61bbb3SSatish Balay           else                  *bap  = value;
10392d61bbb3SSatish Balay           goto noinsert1;
10402d61bbb3SSatish Balay         }
10412d61bbb3SSatish Balay       }
10422d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
10432d61bbb3SSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
10442d61bbb3SSatish Balay       if (nrow >= rmax) {
10452d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
10462d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
10473f1db9ecSBarry Smith         MatScalar *new_a;
10482d61bbb3SSatish Balay 
10492d61bbb3SSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
10502d61bbb3SSatish Balay 
10512d61bbb3SSatish Balay         /* Malloc new storage space */
10523f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
10533f1db9ecSBarry Smith         new_a   = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a);
10542d61bbb3SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
10552d61bbb3SSatish Balay         new_i   = new_j + new_nz;
10562d61bbb3SSatish Balay 
10572d61bbb3SSatish Balay         /* copy over old data into new slots */
10582d61bbb3SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
10592d61bbb3SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1060549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
10612d61bbb3SSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1062549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1063549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1064549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1065549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
10662d61bbb3SSatish Balay         /* free up old matrix storage */
1067606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
1068606d414cSSatish Balay         if (!a->singlemalloc) {
1069606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
1070606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
1071606d414cSSatish Balay         }
10722d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1073c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
10742d61bbb3SSatish Balay 
10752d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10762d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
10773f1db9ecSBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
10782d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10792d61bbb3SSatish Balay         a->reallocs++;
10802d61bbb3SSatish Balay         a->nz++;
10812d61bbb3SSatish Balay       }
10822d61bbb3SSatish Balay       N = nrow++ - 1;
10832d61bbb3SSatish Balay       /* shift up all the later entries in this row */
10842d61bbb3SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
10852d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
1086549d3d68SSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
10872d61bbb3SSatish Balay       }
1088549d3d68SSatish Balay       if (N>=i) {
1089549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1090549d3d68SSatish Balay       }
10912d61bbb3SSatish Balay       rp[i]                      = bcol;
10922d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
10932d61bbb3SSatish Balay       noinsert1:;
10942d61bbb3SSatish Balay       low = i;
10952d61bbb3SSatish Balay     }
10962d61bbb3SSatish Balay     ailen[brow] = nrow;
10972d61bbb3SSatish Balay   }
10982d61bbb3SSatish Balay   PetscFunctionReturn(0);
10992d61bbb3SSatish Balay }
11002d61bbb3SSatish Balay 
11012d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
11022d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
11032d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
11047b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
11057b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
11062d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec);
11072d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
11082d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat);
11092d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *);
11102d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*);
11112d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec);
11122d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
11132d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
11142d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat);
11152d61bbb3SSatish Balay 
11162d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
11172d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
11182d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
11192d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
11202d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
11212d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
112215091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
11232d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
1124f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_7(Mat,Vec,Vec);
1125f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_6(Mat,Vec,Vec);
1126f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_5(Mat,Vec,Vec);
1127f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_4(Mat,Vec,Vec);
1128f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_3(Mat,Vec,Vec);
1129f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_2(Mat,Vec,Vec);
1130f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_1(Mat,Vec,Vec);
11312d61bbb3SSatish Balay 
11322d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
11332d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
11342d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
11352d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
11362d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
11372d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
113815091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
11392d61bbb3SSatish Balay 
11402d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
11412d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
11422d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
11432d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
11442d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
114515091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec);
11462d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
11472d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
11482d61bbb3SSatish Balay 
11492d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
11502d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
11512d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
11522d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
11532d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
115415091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
11552d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
11562d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
11572d61bbb3SSatish Balay 
11582d61bbb3SSatish Balay #undef __FUNC__
11592d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ"
11605ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
11612d61bbb3SSatish Balay {
11622d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
11632d61bbb3SSatish Balay   Mat         outA;
11642d61bbb3SSatish Balay   int         ierr;
1165667159a5SBarry Smith   PetscTruth  row_identity, col_identity;
11662d61bbb3SSatish Balay 
11672d61bbb3SSatish Balay   PetscFunctionBegin;
1168b6d21a55SBarry Smith   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU");
1169667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1170667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1171667159a5SBarry Smith   if (!row_identity || !col_identity) {
1172b008c796SBarry Smith     SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU");
1173667159a5SBarry Smith   }
11742d61bbb3SSatish Balay 
11752d61bbb3SSatish Balay   outA          = inA;
11762d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
11772d61bbb3SSatish Balay 
11782d61bbb3SSatish Balay   if (!a->diag) {
1179c4992f7dSBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
11802d61bbb3SSatish Balay   }
11812d61bbb3SSatish Balay   /*
118215091d37SBarry Smith       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
118315091d37SBarry Smith       for ILU(0) factorization with natural ordering
11842d61bbb3SSatish Balay   */
118515091d37SBarry Smith   switch (a->bs) {
1186f1af5d2fSBarry Smith   case 1:
1187f1af5d2fSBarry Smith     inA->ops->solvetrans      = MatSolveTrans_SeqBAIJ_2_NaturalOrdering;
1188f1af5d2fSBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
118915091d37SBarry Smith   case 2:
119015091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
119115091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_2_NaturalOrdering;
1192f1af5d2fSBarry Smith     inA->ops->solvetrans      = MatSolveTrans_SeqBAIJ_2_NaturalOrdering;
119315091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
119415091d37SBarry Smith     break;
119515091d37SBarry Smith   case 3:
119615091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
119715091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_3_NaturalOrdering;
1198f1af5d2fSBarry Smith     inA->ops->solvetrans      = MatSolveTrans_SeqBAIJ_3_NaturalOrdering;
119915091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
120015091d37SBarry Smith     break;
120115091d37SBarry Smith   case 4:
1202667159a5SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
1203f830108cSBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
1204f1af5d2fSBarry Smith     inA->ops->solvetrans      = MatSolveTrans_SeqBAIJ_4_NaturalOrdering;
120515091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
120615091d37SBarry Smith     break;
120715091d37SBarry Smith   case 5:
1208667159a5SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
1209667159a5SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_5_NaturalOrdering;
1210f1af5d2fSBarry Smith     inA->ops->solvetrans      = MatSolveTrans_SeqBAIJ_5_NaturalOrdering;
121115091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
121215091d37SBarry Smith     break;
121315091d37SBarry Smith   case 6:
121415091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
121515091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_6_NaturalOrdering;
1216f1af5d2fSBarry Smith     inA->ops->solvetrans      = MatSolveTrans_SeqBAIJ_6_NaturalOrdering;
121715091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
121815091d37SBarry Smith     break;
121915091d37SBarry Smith   case 7:
122015091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
1221f1af5d2fSBarry Smith     inA->ops->solvetrans      = MatSolveTrans_SeqBAIJ_7_NaturalOrdering;
122215091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_7_NaturalOrdering;
122315091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
122415091d37SBarry Smith     break;
1225c38d4ed2SBarry Smith   default:
1226c38d4ed2SBarry Smith     a->row        = row;
1227c38d4ed2SBarry Smith     a->col        = col;
1228c38d4ed2SBarry Smith     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1229c38d4ed2SBarry Smith     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1230c38d4ed2SBarry Smith 
1231c38d4ed2SBarry Smith     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1232c38d4ed2SBarry Smith     ierr = ISInvertPermutation(col,&(a->icol));CHKERRQ(ierr);
1233c38d4ed2SBarry Smith     PLogObjectParent(inA,a->icol);
1234c38d4ed2SBarry Smith 
1235c38d4ed2SBarry Smith     if (!a->solve_work) {
1236c38d4ed2SBarry Smith       a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
1237c38d4ed2SBarry Smith       PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
1238c38d4ed2SBarry Smith     }
12392d61bbb3SSatish Balay   }
12402d61bbb3SSatish Balay 
1241667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1242667159a5SBarry Smith 
12432d61bbb3SSatish Balay   PetscFunctionReturn(0);
12442d61bbb3SSatish Balay }
12452d61bbb3SSatish Balay #undef __FUNC__
1246d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ"
1247ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1248ba4ca20aSSatish Balay {
1249c38d4ed2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1250ba4ca20aSSatish Balay   MPI_Comm          comm = A->comm;
1251d132466eSBarry Smith   int               ierr;
1252ba4ca20aSSatish Balay 
12533a40ed3dSBarry Smith   PetscFunctionBegin;
1254c38d4ed2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1255d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1256d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
12573a40ed3dSBarry Smith   PetscFunctionReturn(0);
1258ba4ca20aSSatish Balay }
1259d9b7c43dSSatish Balay 
1260fb2e594dSBarry Smith EXTERN_C_BEGIN
126127a8da17SBarry Smith #undef __FUNC__
126227a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
126327a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
126427a8da17SBarry Smith {
126527a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
126627a8da17SBarry Smith   int         i,nz,n;
126727a8da17SBarry Smith 
126827a8da17SBarry Smith   PetscFunctionBegin;
126927a8da17SBarry Smith   nz = baij->maxnz;
127027a8da17SBarry Smith   n  = baij->n;
127127a8da17SBarry Smith   for (i=0; i<nz; i++) {
127227a8da17SBarry Smith     baij->j[i] = indices[i];
127327a8da17SBarry Smith   }
127427a8da17SBarry Smith   baij->nz = nz;
127527a8da17SBarry Smith   for ( i=0; i<n; i++ ) {
127627a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
127727a8da17SBarry Smith   }
127827a8da17SBarry Smith 
127927a8da17SBarry Smith   PetscFunctionReturn(0);
128027a8da17SBarry Smith }
1281fb2e594dSBarry Smith EXTERN_C_END
128227a8da17SBarry Smith 
128327a8da17SBarry Smith #undef __FUNC__
128427a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices"
128527a8da17SBarry Smith /*@
128627a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
128727a8da17SBarry Smith        in the matrix.
128827a8da17SBarry Smith 
128927a8da17SBarry Smith   Input Parameters:
129027a8da17SBarry Smith +  mat - the SeqBAIJ matrix
129127a8da17SBarry Smith -  indices - the column indices
129227a8da17SBarry Smith 
129315091d37SBarry Smith   Level: advanced
129415091d37SBarry Smith 
129527a8da17SBarry Smith   Notes:
129627a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
129727a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
129827a8da17SBarry Smith   of the MatSetValues() operation.
129927a8da17SBarry Smith 
130027a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
130127a8da17SBarry Smith   MatCreateSeqBAIJ().
130227a8da17SBarry Smith 
130327a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
130427a8da17SBarry Smith 
130527a8da17SBarry Smith @*/
130627a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
130727a8da17SBarry Smith {
130827a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
130927a8da17SBarry Smith 
131027a8da17SBarry Smith   PetscFunctionBegin;
131127a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1312549d3d68SSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr);
131327a8da17SBarry Smith   if (f) {
131427a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
131527a8da17SBarry Smith   } else {
131627a8da17SBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
131727a8da17SBarry Smith   }
131827a8da17SBarry Smith   PetscFunctionReturn(0);
131927a8da17SBarry Smith }
132027a8da17SBarry Smith 
13212593348eSBarry Smith /* -------------------------------------------------------------------*/
1322cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1323cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1324cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1325cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1326cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
1327cc2dc46cSBarry Smith        MatMultTrans_SeqBAIJ,
1328cc2dc46cSBarry Smith        MatMultTransAdd_SeqBAIJ,
1329cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1330cc2dc46cSBarry Smith        0,
1331cc2dc46cSBarry Smith        0,
1332cc2dc46cSBarry Smith        0,
1333cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1334cc2dc46cSBarry Smith        0,
1335b6490206SBarry Smith        0,
1336f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1337cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1338cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1339cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1340cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1341cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1342b6490206SBarry Smith        0,
1343cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1344cc2dc46cSBarry Smith        0,
1345cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1346cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1347cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1348cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1349cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1350cc2dc46cSBarry Smith        0,
1351cc2dc46cSBarry Smith        0,
1352cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1353cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1354cc2dc46cSBarry Smith        MatGetOwnershipRange_SeqBAIJ,
1355cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1356cc2dc46cSBarry Smith        0,
1357cc2dc46cSBarry Smith        0,
1358cc2dc46cSBarry Smith        0,
13592e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1360cc2dc46cSBarry Smith        0,
1361cc2dc46cSBarry Smith        0,
1362cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1363cc2dc46cSBarry Smith        0,
1364cc2dc46cSBarry Smith        0,
1365cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1366cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1367cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1368cc2dc46cSBarry Smith        0,
1369cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1370cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1371cc2dc46cSBarry Smith        0,
1372cc2dc46cSBarry Smith        0,
1373cc2dc46cSBarry Smith        0,
1374cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
13753b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
137692c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1377cc2dc46cSBarry Smith        0,
1378cc2dc46cSBarry Smith        0,
1379cc2dc46cSBarry Smith        0,
1380cc2dc46cSBarry Smith        0,
1381cc2dc46cSBarry Smith        0,
1382cc2dc46cSBarry Smith        0,
1383d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1384cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1385cc2dc46cSBarry Smith        0,
1386cc2dc46cSBarry Smith        0,
1387cc2dc46cSBarry Smith        MatGetMaps_Petsc};
13882593348eSBarry Smith 
13893e90b805SBarry Smith EXTERN_C_BEGIN
13903e90b805SBarry Smith #undef __FUNC__
13913e90b805SBarry Smith #define __FUNC__ "MatStoreValues_SeqBAIJ"
13923e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
13933e90b805SBarry Smith {
13943e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
13953e90b805SBarry Smith   int         nz = aij->i[aij->m]*aij->bs*aij->bs2;
1396d132466eSBarry Smith   int         ierr;
13973e90b805SBarry Smith 
13983e90b805SBarry Smith   PetscFunctionBegin;
13993e90b805SBarry Smith   if (aij->nonew != 1) {
14003e90b805SBarry Smith     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
14013e90b805SBarry Smith   }
14023e90b805SBarry Smith 
14033e90b805SBarry Smith   /* allocate space for values if not already there */
14043e90b805SBarry Smith   if (!aij->saved_values) {
14053e90b805SBarry Smith     aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values);
14063e90b805SBarry Smith   }
14073e90b805SBarry Smith 
14083e90b805SBarry Smith   /* copy values over */
1409d132466eSBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr);
14103e90b805SBarry Smith   PetscFunctionReturn(0);
14113e90b805SBarry Smith }
14123e90b805SBarry Smith EXTERN_C_END
14133e90b805SBarry Smith 
14143e90b805SBarry Smith EXTERN_C_BEGIN
14153e90b805SBarry Smith #undef __FUNC__
141632e87ba3SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqBAIJ"
141732e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
14183e90b805SBarry Smith {
14193e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1420549d3d68SSatish Balay   int         nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr;
14213e90b805SBarry Smith 
14223e90b805SBarry Smith   PetscFunctionBegin;
14233e90b805SBarry Smith   if (aij->nonew != 1) {
14243e90b805SBarry Smith     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
14253e90b805SBarry Smith   }
14263e90b805SBarry Smith   if (!aij->saved_values) {
14273e90b805SBarry Smith     SETERRQ(1,1,"Must call MatStoreValues(A);first");
14283e90b805SBarry Smith   }
14293e90b805SBarry Smith 
14303e90b805SBarry Smith   /* copy values over */
1431549d3d68SSatish Balay   ierr = PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr);
14323e90b805SBarry Smith   PetscFunctionReturn(0);
14333e90b805SBarry Smith }
14343e90b805SBarry Smith EXTERN_C_END
14353e90b805SBarry Smith 
14365615d1e5SSatish Balay #undef __FUNC__
14375615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ"
14382593348eSBarry Smith /*@C
143944cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
144044cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
144144cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
14427fc3c18eSBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
14432bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
14442593348eSBarry Smith 
1445db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1446db81eaa0SLois Curfman McInnes 
14472593348eSBarry Smith    Input Parameters:
1448db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
1449b6490206SBarry Smith .  bs - size of block
14502593348eSBarry Smith .  m - number of rows
14512593348eSBarry Smith .  n - number of columns
1452b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
14537fc3c18eSBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
14542bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
14552593348eSBarry Smith 
14562593348eSBarry Smith    Output Parameter:
14572593348eSBarry Smith .  A - the matrix
14582593348eSBarry Smith 
14590513a670SBarry Smith    Options Database Keys:
1460db81eaa0SLois Curfman McInnes .   -mat_no_unroll - uses code that does not unroll the loops in the
1461db81eaa0SLois Curfman McInnes                      block calculations (much slower)
1462db81eaa0SLois Curfman McInnes .    -mat_block_size - size of the blocks to use
14630513a670SBarry Smith 
146415091d37SBarry Smith    Level: intermediate
146515091d37SBarry Smith 
14662593348eSBarry Smith    Notes:
146744cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
14682593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
146944cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
14702593348eSBarry Smith 
14712593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
14722593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
14732593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
14746da5968aSLois Curfman McInnes    matrices.
14752593348eSBarry Smith 
1476db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
14772593348eSBarry Smith @*/
1478b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
14792593348eSBarry Smith {
14802593348eSBarry Smith   Mat         B;
1481b6490206SBarry Smith   Mat_SeqBAIJ *b;
1482f1af5d2fSBarry Smith   int         i,len,ierr,mbs,nbs,bs2,size;
1483f1af5d2fSBarry Smith   PetscTruth  flg;
14843b2fbd54SBarry Smith 
14853a40ed3dSBarry Smith   PetscFunctionBegin;
1486d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1487a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1488b6490206SBarry Smith 
1489962fb4a0SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1490302e33e3SBarry Smith   mbs  = m/bs;
1491302e33e3SBarry Smith   nbs  = n/bs;
1492302e33e3SBarry Smith   bs2  = bs*bs;
14930513a670SBarry Smith 
14943a40ed3dSBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1495a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
14963a40ed3dSBarry Smith   }
14972593348eSBarry Smith 
1498b73539f3SBarry Smith   if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz);
1499b73539f3SBarry Smith   if (nnz) {
1500faf3f1fcSBarry Smith     for (i=0; i<mbs; i++) {
1501b73539f3SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
1502b73539f3SBarry Smith     }
1503b73539f3SBarry Smith   }
1504b73539f3SBarry Smith 
15052593348eSBarry Smith   *A = 0;
15063f1db9ecSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView);
15072593348eSBarry Smith   PLogObjectCreate(B);
1508b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b);
1509549d3d68SSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1510549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
15111a6a6d98SBarry Smith   ierr    = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
15121a6a6d98SBarry Smith   if (!flg) {
15137fc0212eSBarry Smith     switch (bs) {
15147fc0212eSBarry Smith     case 1:
1515f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1516f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1517f1af5d2fSBarry Smith       B->ops->solvetrans      = MatSolveTrans_SeqBAIJ_1;
1518f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1519f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
15207fc0212eSBarry Smith       break;
15214eeb42bcSBarry Smith     case 2:
1522f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1523f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1524f1af5d2fSBarry Smith       B->ops->solvetrans      = MatSolveTrans_SeqBAIJ_2;
1525f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1526f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
15276d84be18SBarry Smith       break;
1528f327dd97SBarry Smith     case 3:
1529f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1530f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1531f1af5d2fSBarry Smith       B->ops->solvetrans      = MatSolveTrans_SeqBAIJ_3;
1532f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1533f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
15344eeb42bcSBarry Smith       break;
15356d84be18SBarry Smith     case 4:
1536f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1537f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1538f1af5d2fSBarry Smith       B->ops->solvetrans      = MatSolveTrans_SeqBAIJ_4;
1539f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1540f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
15416d84be18SBarry Smith       break;
15426d84be18SBarry Smith     case 5:
1543f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1544f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1545f1af5d2fSBarry Smith       B->ops->solvetrans      = MatSolveTrans_SeqBAIJ_5;
1546f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1547f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
15486d84be18SBarry Smith       break;
154915091d37SBarry Smith     case 6:
155015091d37SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
155115091d37SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_6;
1552f1af5d2fSBarry Smith       B->ops->solvetrans      = MatSolveTrans_SeqBAIJ_6;
155315091d37SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
155415091d37SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
155515091d37SBarry Smith       break;
155648e9ddb2SSatish Balay     case 7:
1557f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1558f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1559f1af5d2fSBarry Smith       B->ops->solvetrans      = MatSolveTrans_SeqBAIJ_7;
1560f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
156148e9ddb2SSatish Balay       break;
15627fc0212eSBarry Smith     }
15631a6a6d98SBarry Smith   }
1564e1311b90SBarry Smith   B->ops->destroy     = MatDestroy_SeqBAIJ;
1565e1311b90SBarry Smith   B->ops->view        = MatView_SeqBAIJ;
15662593348eSBarry Smith   B->factor           = 0;
15672593348eSBarry Smith   B->lupivotthreshold = 1.0;
156890f02eecSBarry Smith   B->mapping          = 0;
15692593348eSBarry Smith   b->row              = 0;
15702593348eSBarry Smith   b->col              = 0;
1571e51c0b9cSSatish Balay   b->icol             = 0;
15722593348eSBarry Smith   b->reallocs         = 0;
15733e90b805SBarry Smith   b->saved_values     = 0;
15742593348eSBarry Smith 
157544cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
157644cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1577a5ae1ecdSBarry Smith 
15787413bad6SBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
15797413bad6SBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1580a5ae1ecdSBarry Smith 
1581b6490206SBarry Smith   b->mbs     = mbs;
1582f2501298SSatish Balay   b->nbs     = nbs;
1583b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(b->imax);
1584c4992f7dSBarry Smith   if (!nnz) {
1585119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
15862593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1587b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1588b6490206SBarry Smith     nz = nz*mbs;
15893a40ed3dSBarry Smith   } else {
15902593348eSBarry Smith     nz = 0;
1591b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
15922593348eSBarry Smith   }
15932593348eSBarry Smith 
15942593348eSBarry Smith   /* allocate the matrix space */
15953f1db9ecSBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int);
15963f1db9ecSBarry Smith   b->a  = (MatScalar *) PetscMalloc( len );CHKPTRQ(b->a);
1597549d3d68SSatish Balay   ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
15987e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
1599549d3d68SSatish Balay   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
16002593348eSBarry Smith   b->i  = b->j + nz;
1601c4992f7dSBarry Smith   b->singlemalloc = PETSC_TRUE;
16022593348eSBarry Smith 
1603de6a44a3SBarry Smith   b->i[0] = 0;
1604b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
16052593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
16062593348eSBarry Smith   }
16072593348eSBarry Smith 
1608b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1609b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1610f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1611b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
16122593348eSBarry Smith 
1613b6490206SBarry Smith   b->bs               = bs;
16149df24120SSatish Balay   b->bs2              = bs2;
1615b6490206SBarry Smith   b->mbs              = mbs;
16162593348eSBarry Smith   b->nz               = 0;
161718e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
1618c4992f7dSBarry Smith   b->sorted           = PETSC_FALSE;
1619c4992f7dSBarry Smith   b->roworiented      = PETSC_TRUE;
16202593348eSBarry Smith   b->nonew            = 0;
16212593348eSBarry Smith   b->diag             = 0;
16222593348eSBarry Smith   b->solve_work       = 0;
1623de6a44a3SBarry Smith   b->mult_work        = 0;
16242593348eSBarry Smith   b->spptr            = 0;
16254e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
1626883fce79SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
16274e220ebcSLois Curfman McInnes 
16282593348eSBarry Smith   *A = B;
16292593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr);
16302593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); }
163127a8da17SBarry Smith 
1632f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
16333e90b805SBarry Smith                                      "MatStoreValues_SeqBAIJ",
16343e90b805SBarry Smith                                      (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1635f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
16363e90b805SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
16373e90b805SBarry Smith                                      (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1638f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
163927a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
164027a8da17SBarry Smith                                      (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
16413a40ed3dSBarry Smith   PetscFunctionReturn(0);
16422593348eSBarry Smith }
16432593348eSBarry Smith 
16445615d1e5SSatish Balay #undef __FUNC__
16452e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ"
16462e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
16472593348eSBarry Smith {
16482593348eSBarry Smith   Mat         C;
1649b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
165027a8da17SBarry Smith   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1651de6a44a3SBarry Smith 
16523a40ed3dSBarry Smith   PetscFunctionBegin;
1653a8c6a408SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
16542593348eSBarry Smith 
16552593348eSBarry Smith   *B = 0;
16563f1db9ecSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView);
16572593348eSBarry Smith   PLogObjectCreate(C);
1658b6490206SBarry Smith   C->data           = (void *) (c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c);
1659549d3d68SSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1660e1311b90SBarry Smith   C->ops->destroy   = MatDestroy_SeqBAIJ;
1661e1311b90SBarry Smith   C->ops->view      = MatView_SeqBAIJ;
16622593348eSBarry Smith   C->factor         = A->factor;
16632593348eSBarry Smith   c->row            = 0;
16642593348eSBarry Smith   c->col            = 0;
1665cac97260SSatish Balay   c->icol           = 0;
166632e87ba3SBarry Smith   c->saved_values   = 0;
1667*f1e2ffcdSBarry Smith   c->keepzeroedrows = a->keepzeroedrows;
16682593348eSBarry Smith   C->assembled      = PETSC_TRUE;
16692593348eSBarry Smith 
167044cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
167144cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
167244cd7ae7SLois Curfman McInnes   C->M          = a->m;
167344cd7ae7SLois Curfman McInnes   C->N          = a->n;
167444cd7ae7SLois Curfman McInnes 
1675b6490206SBarry Smith   c->bs         = a->bs;
16769df24120SSatish Balay   c->bs2        = a->bs2;
1677b6490206SBarry Smith   c->mbs        = a->mbs;
1678b6490206SBarry Smith   c->nbs        = a->nbs;
16792593348eSBarry Smith 
1680b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax);
1681b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen);
1682b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
16832593348eSBarry Smith     c->imax[i] = a->imax[i];
16842593348eSBarry Smith     c->ilen[i] = a->ilen[i];
16852593348eSBarry Smith   }
16862593348eSBarry Smith 
16872593348eSBarry Smith   /* allocate the matrix space */
1688c4992f7dSBarry Smith   c->singlemalloc = PETSC_TRUE;
16893f1db9ecSBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
16903f1db9ecSBarry Smith   c->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(c->a);
16917e67e3f9SSatish Balay   c->j = (int *) (c->a + nz*bs2);
1692de6a44a3SBarry Smith   c->i = c->j + nz;
1693549d3d68SSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1694b6490206SBarry Smith   if (mbs > 0) {
1695549d3d68SSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
16962e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1697549d3d68SSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
16982e8a6d31SBarry Smith     } else {
1699549d3d68SSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
17002593348eSBarry Smith     }
17012593348eSBarry Smith   }
17022593348eSBarry Smith 
1703f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
17042593348eSBarry Smith   c->sorted      = a->sorted;
17052593348eSBarry Smith   c->roworiented = a->roworiented;
17062593348eSBarry Smith   c->nonew       = a->nonew;
17072593348eSBarry Smith 
17082593348eSBarry Smith   if (a->diag) {
1709b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(c->diag);
1710b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1711b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
17122593348eSBarry Smith       c->diag[i] = a->diag[i];
17132593348eSBarry Smith     }
171498305bb5SBarry Smith   } else c->diag        = 0;
17152593348eSBarry Smith   c->nz                 = a->nz;
17162593348eSBarry Smith   c->maxnz              = a->maxnz;
17172593348eSBarry Smith   c->solve_work         = 0;
17182593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
17197fc0212eSBarry Smith   c->mult_work          = 0;
17202593348eSBarry Smith   *B = C;
17217b2a1423SBarry Smith   ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
17223a40ed3dSBarry Smith   PetscFunctionReturn(0);
17232593348eSBarry Smith }
17242593348eSBarry Smith 
17255615d1e5SSatish Balay #undef __FUNC__
17265615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ"
172719bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
17282593348eSBarry Smith {
1729b6490206SBarry Smith   Mat_SeqBAIJ  *a;
17302593348eSBarry Smith   Mat          B;
1731f1af5d2fSBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1732b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
173335aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1734a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1735b6490206SBarry Smith   Scalar       *aa;
173619bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
17372593348eSBarry Smith 
17383a40ed3dSBarry Smith   PetscFunctionBegin;
1739f1af5d2fSBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1740de6a44a3SBarry Smith   bs2  = bs*bs;
1741b6490206SBarry Smith 
1742d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1743cc0fae11SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
174490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
17450752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1746a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
17472593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
17482593348eSBarry Smith 
1749d64ed03dSBarry Smith   if (header[3] < 0) {
1750a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1751d64ed03dSBarry Smith   }
1752d64ed03dSBarry Smith 
1753a8c6a408SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
175435aab85fSBarry Smith 
175535aab85fSBarry Smith   /*
175635aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
175735aab85fSBarry Smith     divisible by the blocksize
175835aab85fSBarry Smith   */
1759b6490206SBarry Smith   mbs        = M/bs;
176035aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
176135aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
176235aab85fSBarry Smith   else                  mbs++;
176335aab85fSBarry Smith   if (extra_rows) {
1764537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
176535aab85fSBarry Smith   }
1766b6490206SBarry Smith 
17672593348eSBarry Smith   /* read in row lengths */
176835aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
17690752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
177035aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
17712593348eSBarry Smith 
1772b6490206SBarry Smith   /* read in column indices */
177335aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) );CHKPTRQ(jj);
17740752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
177535aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1776b6490206SBarry Smith 
1777b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1778b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1779549d3d68SSatish Balay   ierr        = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
178035aab85fSBarry Smith   mask        = (int *) PetscMalloc( 2*mbs*sizeof(int) );CHKPTRQ(mask);
1781549d3d68SSatish Balay   ierr        = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
178235aab85fSBarry Smith   masked      = mask + mbs;
1783b6490206SBarry Smith   rowcount    = 0; nzcount = 0;
1784b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
178535aab85fSBarry Smith     nmask = 0;
1786b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1787b6490206SBarry Smith       kmax = rowlengths[rowcount];
1788b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
178935aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
179035aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1791b6490206SBarry Smith       }
1792b6490206SBarry Smith       rowcount++;
1793b6490206SBarry Smith     }
179435aab85fSBarry Smith     browlengths[i] += nmask;
179535aab85fSBarry Smith     /* zero out the mask elements we set */
179635aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1797b6490206SBarry Smith   }
1798b6490206SBarry Smith 
17992593348eSBarry Smith   /* create our matrix */
18003eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
18012593348eSBarry Smith   B = *A;
1802b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
18032593348eSBarry Smith 
18042593348eSBarry Smith   /* set matrix "i" values */
1805de6a44a3SBarry Smith   a->i[0] = 0;
1806b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1807b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1808b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
18092593348eSBarry Smith   }
1810b6490206SBarry Smith   a->nz         = 0;
1811b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
18122593348eSBarry Smith 
1813b6490206SBarry Smith   /* read in nonzero values */
181435aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
18150752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
181635aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1817b6490206SBarry Smith 
1818b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1819b6490206SBarry Smith   nzcount = 0; jcount = 0;
1820b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1821b6490206SBarry Smith     nzcountb = nzcount;
182235aab85fSBarry Smith     nmask    = 0;
1823b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1824b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1825b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
182635aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
182735aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1828b6490206SBarry Smith       }
1829b6490206SBarry Smith       rowcount++;
1830b6490206SBarry Smith     }
1831de6a44a3SBarry Smith     /* sort the masked values */
1832433994e6SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1833de6a44a3SBarry Smith 
1834b6490206SBarry Smith     /* set "j" values into matrix */
1835b6490206SBarry Smith     maskcount = 1;
183635aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
183735aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1838de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1839b6490206SBarry Smith     }
1840b6490206SBarry Smith     /* set "a" values into matrix */
1841de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1842b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1843b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1844b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1845de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1846de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1847de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1848de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1849b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1850b6490206SBarry Smith       }
1851b6490206SBarry Smith     }
185235aab85fSBarry Smith     /* zero out the mask elements we set */
185335aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1854b6490206SBarry Smith   }
1855a8c6a408SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1856b6490206SBarry Smith 
1857606d414cSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1858606d414cSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1859606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
1860606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
1861606d414cSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
1862b6490206SBarry Smith 
1863b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1864de6a44a3SBarry Smith 
18659c01be13SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
18663a40ed3dSBarry Smith   PetscFunctionReturn(0);
18672593348eSBarry Smith }
18682593348eSBarry Smith 
18692593348eSBarry Smith 
18702593348eSBarry Smith 
1871