xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 1e4dd532b3d0fb4d3279bdf90c27fdefdc7026a5)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*1e4dd532SSatish Balay static char vcid[] = "$Id: baij.c,v 1.128 1998/03/26 22:11:13 balay Exp balay $";
32593348eSBarry Smith #endif
42593348eSBarry Smith 
52593348eSBarry Smith /*
6b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
72593348eSBarry Smith   matrix storage format.
82593348eSBarry Smith */
93369ce9aSBarry Smith 
103369ce9aSBarry Smith #include "pinclude/pviewer.h"
113369ce9aSBarry Smith #include "sys.h"
1270f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
131a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
141a6a6d98SBarry Smith #include "src/inline/spops.h"
152593348eSBarry Smith #include "petsc.h"
163b547af2SSatish Balay 
172d61bbb3SSatish Balay #define CHUNKSIZE  10
18de6a44a3SBarry Smith 
195615d1e5SSatish Balay #undef __FUNC__
205615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ"
21de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A)
22de6a44a3SBarry Smith {
23de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
247fc0212eSBarry Smith   int         i,j, *diag, m = a->mbs;
25de6a44a3SBarry Smith 
263a40ed3dSBarry Smith   PetscFunctionBegin;
27de6a44a3SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
28de6a44a3SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
297fc0212eSBarry Smith   for ( i=0; i<m; i++ ) {
30de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
31de6a44a3SBarry Smith       if (a->j[j] == i) {
32de6a44a3SBarry Smith         diag[i] = j;
33de6a44a3SBarry Smith         break;
34de6a44a3SBarry Smith       }
35de6a44a3SBarry Smith     }
36de6a44a3SBarry Smith   }
37de6a44a3SBarry Smith   a->diag = diag;
383a40ed3dSBarry Smith   PetscFunctionReturn(0);
39de6a44a3SBarry Smith }
402593348eSBarry Smith 
412593348eSBarry Smith 
423b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
433b2fbd54SBarry Smith 
445615d1e5SSatish Balay #undef __FUNC__
455615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ"
463b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
473b2fbd54SBarry Smith                             PetscTruth *done)
483b2fbd54SBarry Smith {
493b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
503b2fbd54SBarry Smith   int         ierr, n = a->mbs,i;
513b2fbd54SBarry Smith 
523a40ed3dSBarry Smith   PetscFunctionBegin;
533b2fbd54SBarry Smith   *nn = n;
543a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
553b2fbd54SBarry Smith   if (symmetric) {
563b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
573b2fbd54SBarry Smith   } else if (oshift == 1) {
583b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
593b2fbd54SBarry Smith     int nz = a->i[n] + 1;
603b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]++;
613b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]++;
623b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
633b2fbd54SBarry Smith   } else {
643b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
653b2fbd54SBarry Smith   }
663b2fbd54SBarry Smith 
673a40ed3dSBarry Smith   PetscFunctionReturn(0);
683b2fbd54SBarry Smith }
693b2fbd54SBarry Smith 
705615d1e5SSatish Balay #undef __FUNC__
71d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ"
723b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
733b2fbd54SBarry Smith                                 PetscTruth *done)
743b2fbd54SBarry Smith {
753b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
763b2fbd54SBarry Smith   int         i,n = a->mbs;
773b2fbd54SBarry Smith 
783a40ed3dSBarry Smith   PetscFunctionBegin;
793a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
803b2fbd54SBarry Smith   if (symmetric) {
813b2fbd54SBarry Smith     PetscFree(*ia);
823b2fbd54SBarry Smith     PetscFree(*ja);
83af5da2bfSBarry Smith   } else if (oshift == 1) {
843b2fbd54SBarry Smith     int nz = a->i[n];
853b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]--;
863b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]--;
873b2fbd54SBarry Smith   }
883a40ed3dSBarry Smith   PetscFunctionReturn(0);
893b2fbd54SBarry Smith }
903b2fbd54SBarry Smith 
912d61bbb3SSatish Balay #undef __FUNC__
922d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ"
932d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs)
942d61bbb3SSatish Balay {
952d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data;
962d61bbb3SSatish Balay 
972d61bbb3SSatish Balay   PetscFunctionBegin;
982d61bbb3SSatish Balay   *bs = baij->bs;
992d61bbb3SSatish Balay   PetscFunctionReturn(0);
1002d61bbb3SSatish Balay }
1012d61bbb3SSatish Balay 
1022d61bbb3SSatish Balay 
1032d61bbb3SSatish Balay #undef __FUNC__
1042d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ"
1052d61bbb3SSatish Balay int MatDestroy_SeqBAIJ(PetscObject obj)
1062d61bbb3SSatish Balay {
1072d61bbb3SSatish Balay   Mat         A  = (Mat) obj;
1082d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
109e51c0b9cSSatish Balay   int         ierr;
1102d61bbb3SSatish Balay 
1112d61bbb3SSatish Balay #if defined(USE_PETSC_LOG)
1122d61bbb3SSatish Balay   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
1132d61bbb3SSatish Balay #endif
1142d61bbb3SSatish Balay   PetscFree(a->a);
1152d61bbb3SSatish Balay   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
1162d61bbb3SSatish Balay   if (a->diag) PetscFree(a->diag);
1172d61bbb3SSatish Balay   if (a->ilen) PetscFree(a->ilen);
1182d61bbb3SSatish Balay   if (a->imax) PetscFree(a->imax);
1192d61bbb3SSatish Balay   if (a->solve_work) PetscFree(a->solve_work);
1202d61bbb3SSatish Balay   if (a->mult_work) PetscFree(a->mult_work);
121e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
1222d61bbb3SSatish Balay   PetscFree(a);
1232d61bbb3SSatish Balay   PLogObjectDestroy(A);
1242d61bbb3SSatish Balay   PetscHeaderDestroy(A);
1252d61bbb3SSatish Balay   PetscFunctionReturn(0);
1262d61bbb3SSatish Balay }
1272d61bbb3SSatish Balay 
1282d61bbb3SSatish Balay #undef __FUNC__
1292d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ"
1302d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1312d61bbb3SSatish Balay {
1322d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1332d61bbb3SSatish Balay 
1342d61bbb3SSatish Balay   PetscFunctionBegin;
1352d61bbb3SSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
1362d61bbb3SSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
1372d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
1382d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
1392d61bbb3SSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
1402d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
1412d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
1422d61bbb3SSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
1432d61bbb3SSatish Balay   else if (op == MAT_ROWS_SORTED ||
1442d61bbb3SSatish Balay            op == MAT_ROWS_UNSORTED ||
1452d61bbb3SSatish Balay            op == MAT_SYMMETRIC ||
1462d61bbb3SSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
1472d61bbb3SSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
148b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
149b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE) {
150981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
1512d61bbb3SSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
1522d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
1532d61bbb3SSatish Balay   } else {
1542d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
1552d61bbb3SSatish Balay   }
1562d61bbb3SSatish Balay   PetscFunctionReturn(0);
1572d61bbb3SSatish Balay }
1582d61bbb3SSatish Balay 
1592d61bbb3SSatish Balay 
1602d61bbb3SSatish Balay #undef __FUNC__
1612d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ"
1622d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
1632d61bbb3SSatish Balay {
1642d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1652d61bbb3SSatish Balay 
1662d61bbb3SSatish Balay   PetscFunctionBegin;
1672d61bbb3SSatish Balay   if (m) *m = a->m;
1682d61bbb3SSatish Balay   if (n) *n = a->n;
1692d61bbb3SSatish Balay   PetscFunctionReturn(0);
1702d61bbb3SSatish Balay }
1712d61bbb3SSatish Balay 
1722d61bbb3SSatish Balay #undef __FUNC__
1732d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ"
1742d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
1752d61bbb3SSatish Balay {
1762d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1772d61bbb3SSatish Balay 
1782d61bbb3SSatish Balay   PetscFunctionBegin;
1792d61bbb3SSatish Balay   *m = 0; *n = a->m;
1802d61bbb3SSatish Balay   PetscFunctionReturn(0);
1812d61bbb3SSatish Balay }
1822d61bbb3SSatish Balay 
1832d61bbb3SSatish Balay #undef __FUNC__
1842d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ"
1852d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
1862d61bbb3SSatish Balay {
1872d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1882d61bbb3SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
1892d61bbb3SSatish Balay   Scalar      *aa,*v_i,*aa_i;
1902d61bbb3SSatish Balay 
1912d61bbb3SSatish Balay   PetscFunctionBegin;
1922d61bbb3SSatish Balay   bs  = a->bs;
1932d61bbb3SSatish Balay   ai  = a->i;
1942d61bbb3SSatish Balay   aj  = a->j;
1952d61bbb3SSatish Balay   aa  = a->a;
1962d61bbb3SSatish Balay   bs2 = a->bs2;
1972d61bbb3SSatish Balay 
1982d61bbb3SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
1992d61bbb3SSatish Balay 
2002d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2012d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2022d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2032d61bbb3SSatish Balay   *nz = bs*M;
2042d61bbb3SSatish Balay 
2052d61bbb3SSatish Balay   if (v) {
2062d61bbb3SSatish Balay     *v = 0;
2072d61bbb3SSatish Balay     if (*nz) {
2082d61bbb3SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
2092d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2102d61bbb3SSatish Balay         v_i  = *v + i*bs;
2112d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2122d61bbb3SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
2132d61bbb3SSatish Balay       }
2142d61bbb3SSatish Balay     }
2152d61bbb3SSatish Balay   }
2162d61bbb3SSatish Balay 
2172d61bbb3SSatish Balay   if (idx) {
2182d61bbb3SSatish Balay     *idx = 0;
2192d61bbb3SSatish Balay     if (*nz) {
2202d61bbb3SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
2212d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2222d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2232d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2242d61bbb3SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
2252d61bbb3SSatish Balay       }
2262d61bbb3SSatish Balay     }
2272d61bbb3SSatish Balay   }
2282d61bbb3SSatish Balay   PetscFunctionReturn(0);
2292d61bbb3SSatish Balay }
2302d61bbb3SSatish Balay 
2312d61bbb3SSatish Balay #undef __FUNC__
2322d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ"
2332d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2342d61bbb3SSatish Balay {
2352d61bbb3SSatish Balay   PetscFunctionBegin;
2362d61bbb3SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
2372d61bbb3SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
2382d61bbb3SSatish Balay   PetscFunctionReturn(0);
2392d61bbb3SSatish Balay }
2402d61bbb3SSatish Balay 
2412d61bbb3SSatish Balay #undef __FUNC__
2422d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ"
2432d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
2442d61bbb3SSatish Balay {
2452d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
2462d61bbb3SSatish Balay   Mat         C;
2472d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
2482d61bbb3SSatish Balay   int         *rows,*cols,bs2=a->bs2;
2492d61bbb3SSatish Balay   Scalar      *array=a->a;
2502d61bbb3SSatish Balay 
2512d61bbb3SSatish Balay   PetscFunctionBegin;
2522d61bbb3SSatish Balay   if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
2532d61bbb3SSatish Balay   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
2542d61bbb3SSatish Balay   PetscMemzero(col,(1+nbs)*sizeof(int));
2552d61bbb3SSatish Balay 
2562d61bbb3SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
2572d61bbb3SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
2582d61bbb3SSatish Balay   PetscFree(col);
2592d61bbb3SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
2602d61bbb3SSatish Balay   cols = rows + bs;
2612d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2622d61bbb3SSatish Balay     cols[0] = i*bs;
2632d61bbb3SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
2642d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
2652d61bbb3SSatish Balay     for ( j=0; j<len; j++ ) {
2662d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
2672d61bbb3SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
2682d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
2692d61bbb3SSatish Balay       array += bs2;
2702d61bbb3SSatish Balay     }
2712d61bbb3SSatish Balay   }
2722d61bbb3SSatish Balay   PetscFree(rows);
2732d61bbb3SSatish Balay 
2742d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2752d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2762d61bbb3SSatish Balay 
2772d61bbb3SSatish Balay   if (B != PETSC_NULL) {
2782d61bbb3SSatish Balay     *B = C;
2792d61bbb3SSatish Balay   } else {
280f830108cSBarry Smith     PetscOps       *Abops;
281f830108cSBarry Smith     struct _MatOps *Aops;
282f830108cSBarry Smith 
2832d61bbb3SSatish Balay     /* This isn't really an in-place transpose */
2842d61bbb3SSatish Balay     PetscFree(a->a);
2852d61bbb3SSatish Balay     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
2862d61bbb3SSatish Balay     if (a->diag) PetscFree(a->diag);
2872d61bbb3SSatish Balay     if (a->ilen) PetscFree(a->ilen);
2882d61bbb3SSatish Balay     if (a->imax) PetscFree(a->imax);
2892d61bbb3SSatish Balay     if (a->solve_work) PetscFree(a->solve_work);
2902d61bbb3SSatish Balay     PetscFree(a);
291f830108cSBarry Smith 
292f830108cSBarry Smith     /*
293f830108cSBarry Smith        This is horrible, horrible code. We need to keep the
294f830108cSBarry Smith       A pointers for the bops and ops but copy everything
295f830108cSBarry Smith       else from C.
296f830108cSBarry Smith     */
297f830108cSBarry Smith     Abops = A->bops;
298f830108cSBarry Smith     Aops  = A->ops;
2992d61bbb3SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
300f830108cSBarry Smith     A->bops = Abops;
301f830108cSBarry Smith     A->ops  = Aops;
302f830108cSBarry Smith 
3032d61bbb3SSatish Balay     PetscHeaderDestroy(C);
3042d61bbb3SSatish Balay   }
3052d61bbb3SSatish Balay   PetscFunctionReturn(0);
3062d61bbb3SSatish Balay }
3072d61bbb3SSatish Balay 
3082d61bbb3SSatish Balay 
3092d61bbb3SSatish Balay 
3103b2fbd54SBarry Smith 
3115615d1e5SSatish Balay #undef __FUNC__
312d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary"
313b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
3142593348eSBarry Smith {
315b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3169df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
317b6490206SBarry Smith   Scalar      *aa;
3182593348eSBarry Smith 
3193a40ed3dSBarry Smith   PetscFunctionBegin;
32090ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
3212593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
3222593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3233b2fbd54SBarry Smith 
3242593348eSBarry Smith   col_lens[1] = a->m;
3252593348eSBarry Smith   col_lens[2] = a->n;
3267e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3272593348eSBarry Smith 
3282593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
329b6490206SBarry Smith   count = 0;
330b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
331b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
332b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
333b6490206SBarry Smith     }
3342593348eSBarry Smith   }
3350752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3362593348eSBarry Smith   PetscFree(col_lens);
3372593348eSBarry Smith 
3382593348eSBarry Smith   /* store column indices (zero start index) */
33966963ce1SSatish Balay   jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj);
340b6490206SBarry Smith   count = 0;
341b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
342b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
343b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
344b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
345b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3462593348eSBarry Smith         }
3472593348eSBarry Smith       }
348b6490206SBarry Smith     }
349b6490206SBarry Smith   }
3500752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr);
351b6490206SBarry Smith   PetscFree(jj);
3522593348eSBarry Smith 
3532593348eSBarry Smith   /* store nonzero values */
35466963ce1SSatish Balay   aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa);
355b6490206SBarry Smith   count = 0;
356b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
357b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
358b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
359b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
3607e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
361b6490206SBarry Smith         }
362b6490206SBarry Smith       }
363b6490206SBarry Smith     }
364b6490206SBarry Smith   }
3650752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
366b6490206SBarry Smith   PetscFree(aa);
3673a40ed3dSBarry Smith   PetscFunctionReturn(0);
3682593348eSBarry Smith }
3692593348eSBarry Smith 
3705615d1e5SSatish Balay #undef __FUNC__
371d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII"
372b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
3732593348eSBarry Smith {
374b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3759df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
3762593348eSBarry Smith   FILE        *fd;
3772593348eSBarry Smith   char        *outputname;
3782593348eSBarry Smith 
3793a40ed3dSBarry Smith   PetscFunctionBegin;
38090ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
3812593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
38290ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
383639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
3847ddc982cSLois Curfman McInnes     fprintf(fd,"  block size is %d\n",bs);
3852593348eSBarry Smith   }
386639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
387a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
3882593348eSBarry Smith   }
389639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
39044cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
39144cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
39244cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
39344cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
39444cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
3953a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
396766eeae4SLois Curfman McInnes           if (imag(a->a[bs2*k + l*bs + j]) > 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0)
39744cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
39844cd7ae7SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
399766eeae4SLois Curfman McInnes           else if (imag(a->a[bs2*k + l*bs + j]) < 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0)
400766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
401766eeae4SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j]));
40244cd7ae7SLois Curfman McInnes           else if (real(a->a[bs2*k + l*bs + j]) != 0.0)
40344cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
40444cd7ae7SLois Curfman McInnes #else
40544cd7ae7SLois Curfman McInnes           if (a->a[bs2*k + l*bs + j] != 0.0)
40644cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
40744cd7ae7SLois Curfman McInnes #endif
40844cd7ae7SLois Curfman McInnes           }
40944cd7ae7SLois Curfman McInnes         }
41044cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
41144cd7ae7SLois Curfman McInnes       }
41244cd7ae7SLois Curfman McInnes     }
41344cd7ae7SLois Curfman McInnes   }
4142593348eSBarry Smith   else {
415b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
416b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
417b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
418b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
419b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
4203a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
421766eeae4SLois Curfman McInnes           if (imag(a->a[bs2*k + l*bs + j]) > 0.0) {
42288685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
4237e67e3f9SSatish Balay               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
42488685aaeSLois Curfman McInnes           }
425766eeae4SLois Curfman McInnes           else if (imag(a->a[bs2*k + l*bs + j]) < 0.0) {
426766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
427766eeae4SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j]));
428766eeae4SLois Curfman McInnes           }
42988685aaeSLois Curfman McInnes           else {
4307e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
43188685aaeSLois Curfman McInnes           }
43288685aaeSLois Curfman McInnes #else
4337e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
43488685aaeSLois Curfman McInnes #endif
4352593348eSBarry Smith           }
4362593348eSBarry Smith         }
4372593348eSBarry Smith         fprintf(fd,"\n");
4382593348eSBarry Smith       }
4392593348eSBarry Smith     }
440b6490206SBarry Smith   }
4412593348eSBarry Smith   fflush(fd);
4423a40ed3dSBarry Smith   PetscFunctionReturn(0);
4432593348eSBarry Smith }
4442593348eSBarry Smith 
4455615d1e5SSatish Balay #undef __FUNC__
446d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw"
4473270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
4483270192aSSatish Balay {
4493270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
4503270192aSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2;
4513270192aSSatish Balay   double       xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
4523270192aSSatish Balay   Scalar       *aa;
4533270192aSSatish Balay   Draw         draw;
4543270192aSSatish Balay   DrawButton   button;
4553270192aSSatish Balay   PetscTruth   isnull;
4563270192aSSatish Balay 
4573a40ed3dSBarry Smith   PetscFunctionBegin;
4583a40ed3dSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr);
4593a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
4603270192aSSatish Balay 
4613270192aSSatish Balay   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
4623270192aSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
4633270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
4643270192aSSatish Balay   /* loop over matrix elements drawing boxes */
4653270192aSSatish Balay   color = DRAW_BLUE;
4663270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4673270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4683270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4693270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4703270192aSSatish Balay       aa = a->a + j*bs2;
4713270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4723270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4733270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
474b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4753270192aSSatish Balay         }
4763270192aSSatish Balay       }
4773270192aSSatish Balay     }
4783270192aSSatish Balay   }
4793270192aSSatish Balay   color = DRAW_CYAN;
4803270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4813270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4823270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4833270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4843270192aSSatish Balay       aa = a->a + j*bs2;
4853270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4863270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4873270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
488b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4893270192aSSatish Balay         }
4903270192aSSatish Balay       }
4913270192aSSatish Balay     }
4923270192aSSatish Balay   }
4933270192aSSatish Balay 
4943270192aSSatish Balay   color = DRAW_RED;
4953270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4963270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4973270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4983270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4993270192aSSatish Balay       aa = a->a + j*bs2;
5003270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5013270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5023270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
503b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5043270192aSSatish Balay         }
5053270192aSSatish Balay       }
5063270192aSSatish Balay     }
5073270192aSSatish Balay   }
5083270192aSSatish Balay 
50955843e3eSBarry Smith   DrawSynchronizedFlush(draw);
5103270192aSSatish Balay   DrawGetPause(draw,&pause);
5113a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
5123270192aSSatish Balay 
5133270192aSSatish Balay   /* allow the matrix to zoom or shrink */
51455843e3eSBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
5153270192aSSatish Balay   while (button != BUTTON_RIGHT) {
51655843e3eSBarry Smith     DrawSynchronizedClear(draw);
5173270192aSSatish Balay     if (button == BUTTON_LEFT) scale = .5;
5183270192aSSatish Balay     else if (button == BUTTON_CENTER) scale = 2.;
5193270192aSSatish Balay     xl = scale*(xl + w - xc) + xc - w*scale;
5203270192aSSatish Balay     xr = scale*(xr - w - xc) + xc + w*scale;
5213270192aSSatish Balay     yl = scale*(yl + h - yc) + yc - h*scale;
5223270192aSSatish Balay     yr = scale*(yr - h - yc) + yc + h*scale;
5233270192aSSatish Balay     w *= scale; h *= scale;
5243270192aSSatish Balay     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5253270192aSSatish Balay     color = DRAW_BLUE;
5263270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5273270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5283270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5293270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5303270192aSSatish Balay         aa = a->a + j*bs2;
5313270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5323270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5333270192aSSatish Balay             if (PetscReal(*aa++) >=  0.) continue;
534b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5353270192aSSatish Balay           }
5363270192aSSatish Balay         }
5373270192aSSatish Balay       }
5383270192aSSatish Balay     }
5393270192aSSatish Balay     color = DRAW_CYAN;
5403270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5413270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5423270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5433270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5443270192aSSatish Balay         aa = a->a + j*bs2;
5453270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5463270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5473270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
548b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5493270192aSSatish Balay           }
5503270192aSSatish Balay         }
5513270192aSSatish Balay       }
5523270192aSSatish Balay     }
5533270192aSSatish Balay 
5543270192aSSatish Balay     color = DRAW_RED;
5553270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5563270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5573270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5583270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5593270192aSSatish Balay         aa = a->a + j*bs2;
5603270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5613270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5623270192aSSatish Balay             if (PetscReal(*aa++) <= 0.) continue;
563b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5643270192aSSatish Balay           }
5653270192aSSatish Balay         }
5663270192aSSatish Balay       }
5673270192aSSatish Balay     }
56855843e3eSBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
5693270192aSSatish Balay   }
5703a40ed3dSBarry Smith   PetscFunctionReturn(0);
5713270192aSSatish Balay }
5723270192aSSatish Balay 
5735615d1e5SSatish Balay #undef __FUNC__
574d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ"
5758f6be9afSLois Curfman McInnes int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
5762593348eSBarry Smith {
5772593348eSBarry Smith   Mat         A = (Mat) obj;
57819bcc07fSBarry Smith   ViewerType  vtype;
57919bcc07fSBarry Smith   int         ierr;
5802593348eSBarry Smith 
5813a40ed3dSBarry Smith   PetscFunctionBegin;
58219bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
58319bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
584a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported");
5853a40ed3dSBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
5863a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
5873a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
5883a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
5893a40ed3dSBarry Smith   } else if (vtype == DRAW_VIEWER) {
5903a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
5912593348eSBarry Smith   }
5923a40ed3dSBarry Smith   PetscFunctionReturn(0);
5932593348eSBarry Smith }
594b6490206SBarry Smith 
595cd0e1443SSatish Balay 
5965615d1e5SSatish Balay #undef __FUNC__
5972d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ"
5982d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
599cd0e1443SSatish Balay {
600cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6012d61bbb3SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
6022d61bbb3SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
6032d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6042d61bbb3SSatish Balay   Scalar     *ap, *aa = a->a, zero = 0.0;
605cd0e1443SSatish Balay 
6063a40ed3dSBarry Smith   PetscFunctionBegin;
6072d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
608cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
609a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
610a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
6112d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6122c3acbe9SBarry Smith     nrow = ailen[brow];
6132d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
614a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
615a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6162d61bbb3SSatish Balay       col  = in[l] ;
6172d61bbb3SSatish Balay       bcol = col/bs;
6182d61bbb3SSatish Balay       cidx = col%bs;
6192d61bbb3SSatish Balay       ridx = row%bs;
6202d61bbb3SSatish Balay       high = nrow;
6212d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6222d61bbb3SSatish Balay       while (high-low > 5) {
623cd0e1443SSatish Balay         t = (low+high)/2;
624cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
625cd0e1443SSatish Balay         else             low  = t;
626cd0e1443SSatish Balay       }
627cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
628cd0e1443SSatish Balay         if (rp[i] > bcol) break;
629cd0e1443SSatish Balay         if (rp[i] == bcol) {
6302d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6312d61bbb3SSatish Balay           goto finished;
632cd0e1443SSatish Balay         }
633cd0e1443SSatish Balay       }
6342d61bbb3SSatish Balay       *v++ = zero;
6352d61bbb3SSatish Balay       finished:;
636cd0e1443SSatish Balay     }
637cd0e1443SSatish Balay   }
6383a40ed3dSBarry Smith   PetscFunctionReturn(0);
639cd0e1443SSatish Balay }
640cd0e1443SSatish Balay 
6412d61bbb3SSatish Balay 
6425615d1e5SSatish Balay #undef __FUNC__
64305a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ"
64492c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
64592c4ed94SBarry Smith {
64692c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6478a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
648dd9472c6SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented;
649dd9472c6SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
6500e324ae4SSatish Balay   Scalar      *ap,*value=v,*aa=a->a,*bap;
65192c4ed94SBarry Smith 
6523a40ed3dSBarry Smith   PetscFunctionBegin;
6530e324ae4SSatish Balay   if (roworiented) {
6540e324ae4SSatish Balay     stepval = (n-1)*bs;
6550e324ae4SSatish Balay   } else {
6560e324ae4SSatish Balay     stepval = (m-1)*bs;
6570e324ae4SSatish Balay   }
65892c4ed94SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
65992c4ed94SBarry Smith     row  = im[k];
6603a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
661a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
662a8c6a408SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
66392c4ed94SBarry Smith #endif
66492c4ed94SBarry Smith     rp   = aj + ai[row];
66592c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
66692c4ed94SBarry Smith     rmax = imax[row];
66792c4ed94SBarry Smith     nrow = ailen[row];
66892c4ed94SBarry Smith     low  = 0;
66992c4ed94SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
6703a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
671a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
672a8c6a408SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
67392c4ed94SBarry Smith #endif
67492c4ed94SBarry Smith       col = in[l];
67592c4ed94SBarry Smith       if (roworiented) {
6760e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
6770e324ae4SSatish Balay       } else {
6780e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
67992c4ed94SBarry Smith       }
68092c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
68192c4ed94SBarry Smith       while (high-low > 7) {
68292c4ed94SBarry Smith         t = (low+high)/2;
68392c4ed94SBarry Smith         if (rp[t] > col) high = t;
68492c4ed94SBarry Smith         else             low  = t;
68592c4ed94SBarry Smith       }
68692c4ed94SBarry Smith       for ( i=low; i<high; i++ ) {
68792c4ed94SBarry Smith         if (rp[i] > col) break;
68892c4ed94SBarry Smith         if (rp[i] == col) {
6898a84c255SSatish Balay           bap  = ap +  bs2*i;
6900e324ae4SSatish Balay           if (roworiented) {
6918a84c255SSatish Balay             if (is == ADD_VALUES) {
692dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
693dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
6948a84c255SSatish Balay                   bap[jj] += *value++;
695dd9472c6SBarry Smith                 }
696dd9472c6SBarry Smith               }
6970e324ae4SSatish Balay             } else {
698dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
699dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7000e324ae4SSatish Balay                   bap[jj] = *value++;
7018a84c255SSatish Balay                 }
702dd9472c6SBarry Smith               }
703dd9472c6SBarry Smith             }
7040e324ae4SSatish Balay           } else {
7050e324ae4SSatish Balay             if (is == ADD_VALUES) {
706dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
707dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7080e324ae4SSatish Balay                   *bap++ += *value++;
709dd9472c6SBarry Smith                 }
710dd9472c6SBarry Smith               }
7110e324ae4SSatish Balay             } else {
712dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
713dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7140e324ae4SSatish Balay                   *bap++  = *value++;
7150e324ae4SSatish Balay                 }
7168a84c255SSatish Balay               }
717dd9472c6SBarry Smith             }
718dd9472c6SBarry Smith           }
719f1241b54SBarry Smith           goto noinsert2;
72092c4ed94SBarry Smith         }
72192c4ed94SBarry Smith       }
72289280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
723a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
72492c4ed94SBarry Smith       if (nrow >= rmax) {
72592c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
72692c4ed94SBarry Smith         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
72792c4ed94SBarry Smith         Scalar *new_a;
72892c4ed94SBarry Smith 
729a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
73089280ab3SLois Curfman McInnes 
73192c4ed94SBarry Smith         /* malloc new storage space */
73292c4ed94SBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
73392c4ed94SBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
73492c4ed94SBarry Smith         new_j   = (int *) (new_a + bs2*new_nz);
73592c4ed94SBarry Smith         new_i   = new_j + new_nz;
73692c4ed94SBarry Smith 
73792c4ed94SBarry Smith         /* copy over old data into new slots */
73892c4ed94SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
73992c4ed94SBarry Smith         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
74092c4ed94SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));
74192c4ed94SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow);
742dd9472c6SBarry Smith         PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));
74392c4ed94SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(Scalar));
74492c4ed94SBarry Smith         PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
745dd9472c6SBarry Smith         PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(Scalar));
74692c4ed94SBarry Smith         /* free up old matrix storage */
74792c4ed94SBarry Smith         PetscFree(a->a);
74892c4ed94SBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
74992c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
75092c4ed94SBarry Smith         a->singlemalloc = 1;
75192c4ed94SBarry Smith 
75292c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
75392c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
75492c4ed94SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
75592c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
75692c4ed94SBarry Smith         a->reallocs++;
75792c4ed94SBarry Smith         a->nz++;
75892c4ed94SBarry Smith       }
75992c4ed94SBarry Smith       N = nrow++ - 1;
76092c4ed94SBarry Smith       /* shift up all the later entries in this row */
76192c4ed94SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
76292c4ed94SBarry Smith         rp[ii+1] = rp[ii];
76392c4ed94SBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
76492c4ed94SBarry Smith       }
76592c4ed94SBarry Smith       if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
76692c4ed94SBarry Smith       rp[i] = col;
7678a84c255SSatish Balay       bap   = ap +  bs2*i;
7680e324ae4SSatish Balay       if (roworiented) {
769dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval) {
770dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs ) {
7710e324ae4SSatish Balay             bap[jj] = *value++;
772dd9472c6SBarry Smith           }
773dd9472c6SBarry Smith         }
7740e324ae4SSatish Balay       } else {
775dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval ) {
776dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++ ) {
7770e324ae4SSatish Balay             *bap++  = *value++;
7780e324ae4SSatish Balay           }
779dd9472c6SBarry Smith         }
780dd9472c6SBarry Smith       }
781f1241b54SBarry Smith       noinsert2:;
78292c4ed94SBarry Smith       low = i;
78392c4ed94SBarry Smith     }
78492c4ed94SBarry Smith     ailen[row] = nrow;
78592c4ed94SBarry Smith   }
7863a40ed3dSBarry Smith   PetscFunctionReturn(0);
78792c4ed94SBarry Smith }
78892c4ed94SBarry Smith 
789f2501298SSatish Balay 
7905615d1e5SSatish Balay #undef __FUNC__
7915615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ"
7928f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
793584200bdSSatish Balay {
794584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
795584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
796a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
79743ee02c3SBarry Smith   int        mbs = a->mbs, bs2 = a->bs2,rmax = 0;
798584200bdSSatish Balay   Scalar     *aa = a->a, *ap;
799584200bdSSatish Balay 
8003a40ed3dSBarry Smith   PetscFunctionBegin;
8013a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
802584200bdSSatish Balay 
80343ee02c3SBarry Smith   if (m) rmax = ailen[0];
804584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
805584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
806584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
807d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
808584200bdSSatish Balay     if (fshift) {
809a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
810584200bdSSatish Balay       N = ailen[i];
811584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
812584200bdSSatish Balay         ip[j-fshift] = ip[j];
8137e67e3f9SSatish Balay         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar));
814584200bdSSatish Balay       }
815584200bdSSatish Balay     }
816584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
817584200bdSSatish Balay   }
818584200bdSSatish Balay   if (mbs) {
819584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
820584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
821584200bdSSatish Balay   }
822584200bdSSatish Balay   /* reset ilen and imax for each row */
823584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
824584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
825584200bdSSatish Balay   }
826a7c10996SSatish Balay   a->nz = ai[mbs];
827584200bdSSatish Balay 
828584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
829584200bdSSatish Balay   if (fshift && a->diag) {
830584200bdSSatish Balay     PetscFree(a->diag);
831584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
832584200bdSSatish Balay     a->diag = 0;
833584200bdSSatish Balay   }
8343d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
835219d9a1aSLois Curfman McInnes            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
8363d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
837584200bdSSatish Balay            a->reallocs);
838d402145bSBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
839e2f3b5e9SSatish Balay   a->reallocs          = 0;
8404e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
8414e220ebcSLois Curfman McInnes 
8423a40ed3dSBarry Smith   PetscFunctionReturn(0);
843584200bdSSatish Balay }
844584200bdSSatish Balay 
8455a838052SSatish Balay 
846d9b7c43dSSatish Balay /* idx should be of length atlease bs */
8475615d1e5SSatish Balay #undef __FUNC__
8485615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block"
849d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg)
850d9b7c43dSSatish Balay {
851d9b7c43dSSatish Balay   int i,row;
8523a40ed3dSBarry Smith 
8533a40ed3dSBarry Smith   PetscFunctionBegin;
854d9b7c43dSSatish Balay   row = idx[0];
8553a40ed3dSBarry Smith   if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
856d9b7c43dSSatish Balay 
857d9b7c43dSSatish Balay   for ( i=1; i<bs; i++ ) {
8583a40ed3dSBarry Smith     if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
859d9b7c43dSSatish Balay   }
860d9b7c43dSSatish Balay   *flg = PETSC_TRUE;
8613a40ed3dSBarry Smith   PetscFunctionReturn(0);
862d9b7c43dSSatish Balay }
863d9b7c43dSSatish Balay 
8645615d1e5SSatish Balay #undef __FUNC__
8655615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ"
8668f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
867d9b7c43dSSatish Balay {
868d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
869d9b7c43dSSatish Balay   IS          is_local;
870d9b7c43dSSatish Balay   int         ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2;
871d9b7c43dSSatish Balay   PetscTruth  flg;
872d9b7c43dSSatish Balay   Scalar      zero = 0.0,*aa;
873d9b7c43dSSatish Balay 
8743a40ed3dSBarry Smith   PetscFunctionBegin;
875d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
876d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
877d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
878537820f0SBarry Smith   ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr);
879d9b7c43dSSatish Balay   ierr = ISSort(is_local); CHKERRQ(ierr);
880d9b7c43dSSatish Balay   ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr);
881d9b7c43dSSatish Balay 
882d9b7c43dSSatish Balay   i = 0;
883d9b7c43dSSatish Balay   while (i < is_n) {
884a8c6a408SBarry Smith     if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
885d9b7c43dSSatish Balay     flg = PETSC_FALSE;
886d9b7c43dSSatish Balay     if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); }
887d9b7c43dSSatish Balay     count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs;
888d9b7c43dSSatish Balay     aa    = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs);
8892d61bbb3SSatish Balay     if (flg) { /* There exists a block of rows to be Zerowed */
890a07cd24cSSatish Balay       if (baij->ilen[rows[i]/bs] > 0) {
8912d61bbb3SSatish Balay         PetscMemzero(aa,count*bs*sizeof(Scalar));
8922d61bbb3SSatish Balay         baij->ilen[rows[i]/bs] = 1;
8932d61bbb3SSatish Balay         baij->j[baij->i[rows[i]/bs]] = rows[i]/bs;
894a07cd24cSSatish Balay       }
8952d61bbb3SSatish Balay       i += bs;
8962d61bbb3SSatish Balay     } else { /* Zero out only the requested row */
897d9b7c43dSSatish Balay       for ( j=0; j<count; j++ ) {
898d9b7c43dSSatish Balay         aa[0] = zero;
899d9b7c43dSSatish Balay         aa+=bs;
900d9b7c43dSSatish Balay       }
901d9b7c43dSSatish Balay       i++;
902d9b7c43dSSatish Balay     }
903d9b7c43dSSatish Balay   }
904d9b7c43dSSatish Balay   if (diag) {
905d9b7c43dSSatish Balay     for ( j=0; j<is_n; j++ ) {
906f830108cSBarry Smith       ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
9072d61bbb3SSatish Balay       /* ierr = MatSetValues(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); */
908d9b7c43dSSatish Balay     }
909d9b7c43dSSatish Balay   }
910d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr);
911d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr);
912d9b7c43dSSatish Balay   ierr = ISDestroy(is_local); CHKERRQ(ierr);
9139a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9143a40ed3dSBarry Smith   PetscFunctionReturn(0);
915d9b7c43dSSatish Balay }
9161c351548SSatish Balay 
9175615d1e5SSatish Balay #undef __FUNC__
9182d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ"
9192d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
9202d61bbb3SSatish Balay {
9212d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9222d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
9232d61bbb3SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
9242d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
9252d61bbb3SSatish Balay   int          ridx,cidx,bs2=a->bs2;
9262d61bbb3SSatish Balay   Scalar      *ap,value,*aa=a->a,*bap;
9272d61bbb3SSatish Balay 
9282d61bbb3SSatish Balay   PetscFunctionBegin;
9292d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
9302d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
9312d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9322d61bbb3SSatish Balay     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
9332d61bbb3SSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
9342d61bbb3SSatish Balay #endif
9352d61bbb3SSatish Balay     rp   = aj + ai[brow];
9362d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
9372d61bbb3SSatish Balay     rmax = imax[brow];
9382d61bbb3SSatish Balay     nrow = ailen[brow];
9392d61bbb3SSatish Balay     low  = 0;
9402d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
9412d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9422d61bbb3SSatish Balay       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
9432d61bbb3SSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
9442d61bbb3SSatish Balay #endif
9452d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
9462d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
9472d61bbb3SSatish Balay       if (roworiented) {
9482d61bbb3SSatish Balay         value = *v++;
9492d61bbb3SSatish Balay       } else {
9502d61bbb3SSatish Balay         value = v[k + l*m];
9512d61bbb3SSatish Balay       }
9522d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
9532d61bbb3SSatish Balay       while (high-low > 7) {
9542d61bbb3SSatish Balay         t = (low+high)/2;
9552d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
9562d61bbb3SSatish Balay         else              low  = t;
9572d61bbb3SSatish Balay       }
9582d61bbb3SSatish Balay       for ( i=low; i<high; i++ ) {
9592d61bbb3SSatish Balay         if (rp[i] > bcol) break;
9602d61bbb3SSatish Balay         if (rp[i] == bcol) {
9612d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
9622d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
9632d61bbb3SSatish Balay           else                  *bap  = value;
9642d61bbb3SSatish Balay           goto noinsert1;
9652d61bbb3SSatish Balay         }
9662d61bbb3SSatish Balay       }
9672d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
9682d61bbb3SSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9692d61bbb3SSatish Balay       if (nrow >= rmax) {
9702d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
9712d61bbb3SSatish Balay         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
9722d61bbb3SSatish Balay         Scalar *new_a;
9732d61bbb3SSatish Balay 
9742d61bbb3SSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9752d61bbb3SSatish Balay 
9762d61bbb3SSatish Balay         /* Malloc new storage space */
9772d61bbb3SSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
9782d61bbb3SSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
9792d61bbb3SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
9802d61bbb3SSatish Balay         new_i   = new_j + new_nz;
9812d61bbb3SSatish Balay 
9822d61bbb3SSatish Balay         /* copy over old data into new slots */
9832d61bbb3SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
9842d61bbb3SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
9852d61bbb3SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
9862d61bbb3SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
9872d61bbb3SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
9882d61bbb3SSatish Balay                                                            len*sizeof(int));
9892d61bbb3SSatish Balay         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar));
9902d61bbb3SSatish Balay         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
9912d61bbb3SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
9922d61bbb3SSatish Balay                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar));
9932d61bbb3SSatish Balay         /* free up old matrix storage */
9942d61bbb3SSatish Balay         PetscFree(a->a);
9952d61bbb3SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
9962d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
9972d61bbb3SSatish Balay         a->singlemalloc = 1;
9982d61bbb3SSatish Balay 
9992d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10002d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
10012d61bbb3SSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
10022d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10032d61bbb3SSatish Balay         a->reallocs++;
10042d61bbb3SSatish Balay         a->nz++;
10052d61bbb3SSatish Balay       }
10062d61bbb3SSatish Balay       N = nrow++ - 1;
10072d61bbb3SSatish Balay       /* shift up all the later entries in this row */
10082d61bbb3SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
10092d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
10102d61bbb3SSatish Balay         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
10112d61bbb3SSatish Balay       }
10122d61bbb3SSatish Balay       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
10132d61bbb3SSatish Balay       rp[i]                      = bcol;
10142d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
10152d61bbb3SSatish Balay       noinsert1:;
10162d61bbb3SSatish Balay       low = i;
10172d61bbb3SSatish Balay     }
10182d61bbb3SSatish Balay     ailen[brow] = nrow;
10192d61bbb3SSatish Balay   }
10202d61bbb3SSatish Balay   PetscFunctionReturn(0);
10212d61bbb3SSatish Balay }
10222d61bbb3SSatish Balay 
10232d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
10242d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
10252d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
1026d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*);
10272d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**);
10282d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec);
10292d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
10302d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat);
10312d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *);
10322d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*);
10332d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec);
10342d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
10352d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
10362d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat);
10372d61bbb3SSatish Balay 
10382d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
10392d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
10402d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
10412d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
10422d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
10432d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
10442d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
10452d61bbb3SSatish Balay 
10462d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
10472d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
10482d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
10492d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
10502d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
10512d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
10522d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
10532d61bbb3SSatish Balay 
10542d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
10552d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
10562d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
10572d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
10582d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
10592d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
10602d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
10612d61bbb3SSatish Balay 
10622d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
10632d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
10642d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
10652d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
10662d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
10672d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
10682d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
10692d61bbb3SSatish Balay 
10702d61bbb3SSatish Balay /*
10712d61bbb3SSatish Balay      note: This can only work for identity for row and col. It would
10722d61bbb3SSatish Balay    be good to check this and otherwise generate an error.
10732d61bbb3SSatish Balay */
10742d61bbb3SSatish Balay #undef __FUNC__
10752d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ"
10762d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
10772d61bbb3SSatish Balay {
10782d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
10792d61bbb3SSatish Balay   Mat         outA;
10802d61bbb3SSatish Balay   int         ierr;
10812d61bbb3SSatish Balay 
10822d61bbb3SSatish Balay   PetscFunctionBegin;
10832d61bbb3SSatish Balay   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
10842d61bbb3SSatish Balay 
10852d61bbb3SSatish Balay   outA          = inA;
10862d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
10872d61bbb3SSatish Balay   a->row        = row;
10882d61bbb3SSatish Balay   a->col        = col;
10892d61bbb3SSatish Balay 
1090e51c0b9cSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1091e51c0b9cSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
1092*1e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1093e51c0b9cSSatish Balay 
10942d61bbb3SSatish Balay   if (!a->solve_work) {
10952d61bbb3SSatish Balay     a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
10962d61bbb3SSatish Balay     PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
10972d61bbb3SSatish Balay   }
10982d61bbb3SSatish Balay 
10992d61bbb3SSatish Balay   if (!a->diag) {
11002d61bbb3SSatish Balay     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
11012d61bbb3SSatish Balay   }
11022d61bbb3SSatish Balay   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
11032d61bbb3SSatish Balay 
11042d61bbb3SSatish Balay   /*
11052d61bbb3SSatish Balay       Blocksize 4 has a special faster solver for ILU(0) factorization
11062d61bbb3SSatish Balay     with natural ordering
11072d61bbb3SSatish Balay   */
11082d61bbb3SSatish Balay   if (a->bs == 4) {
1109f830108cSBarry Smith     inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering;
11102d61bbb3SSatish Balay   }
11112d61bbb3SSatish Balay 
11122d61bbb3SSatish Balay   PetscFunctionReturn(0);
11132d61bbb3SSatish Balay }
11142d61bbb3SSatish Balay #undef __FUNC__
1115d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ"
1116ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1117ba4ca20aSSatish Balay {
1118ba4ca20aSSatish Balay   static int called = 0;
1119ba4ca20aSSatish Balay   MPI_Comm   comm = A->comm;
1120ba4ca20aSSatish Balay 
11213a40ed3dSBarry Smith   PetscFunctionBegin;
11223a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
112376be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");
112476be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");
11253a40ed3dSBarry Smith   PetscFunctionReturn(0);
1126ba4ca20aSSatish Balay }
1127d9b7c43dSSatish Balay 
11282593348eSBarry Smith /* -------------------------------------------------------------------*/
1129cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ,
11309867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
1131f44d3a6dSSatish Balay        MatMult_SeqBAIJ_N,MatMultAdd_SeqBAIJ_N,
1132faf6580fSSatish Balay        MatMultTrans_SeqBAIJ,MatMultTransAdd_SeqBAIJ,
11331a6a6d98SBarry Smith        MatSolve_SeqBAIJ_N,0,
1134b6490206SBarry Smith        0,0,
1135de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
1136b6490206SBarry Smith        0,
1137f2501298SSatish Balay        MatTranspose_SeqBAIJ,
113817e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
11399867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
1140584200bdSSatish Balay        0,MatAssemblyEnd_SeqBAIJ,
1141b6490206SBarry Smith        0,
1142d9b7c43dSSatish Balay        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,MatZeroRows_SeqBAIJ,
11437fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
1144b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
1145de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
1146d402145bSBarry Smith        0,0,
1147b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
1148b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
1149af015674SSatish Balay        MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ,
11507e67e3f9SSatish Balay        MatGetValues_SeqBAIJ,0,
1151ba4ca20aSSatish Balay        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
11523b2fbd54SBarry Smith        0,0,0,MatGetBlockSize_SeqBAIJ,
11533b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
115492c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
115592c4ed94SBarry Smith        0,0,0,0,0,0,
1156d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1157d3825aa8SBarry Smith        MatGetSubMatrix_SeqBAIJ};
11582593348eSBarry Smith 
11595615d1e5SSatish Balay #undef __FUNC__
11605615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ"
11612593348eSBarry Smith /*@C
116244cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
116344cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
116444cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
11652bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
11662bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
11672593348eSBarry Smith 
11682593348eSBarry Smith    Input Parameters:
1169029af93fSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
1170b6490206SBarry Smith .  bs - size of block
11712593348eSBarry Smith .  m - number of rows
11722593348eSBarry Smith .  n - number of columns
1173b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
11742bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of block nonzeros in the various block rows
11752bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
11762593348eSBarry Smith 
11772593348eSBarry Smith    Output Parameter:
11782593348eSBarry Smith .  A - the matrix
11792593348eSBarry Smith 
11800513a670SBarry Smith    Options Database Keys:
11810513a670SBarry Smith $    -mat_no_unroll - uses code that does not unroll the loops in the
11820effe34fSLois Curfman McInnes $                     block calculations (much slower)
11830513a670SBarry Smith $    -mat_block_size - size of the blocks to use
11840513a670SBarry Smith 
11852593348eSBarry Smith    Notes:
118644cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
11872593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
118844cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
11892593348eSBarry Smith 
11902593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
11912593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
11922593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
11936da5968aSLois Curfman McInnes    matrices.
11942593348eSBarry Smith 
119544cd7ae7SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
11962593348eSBarry Smith @*/
1197b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
11982593348eSBarry Smith {
11992593348eSBarry Smith   Mat         B;
1200b6490206SBarry Smith   Mat_SeqBAIJ *b;
12013b2fbd54SBarry Smith   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size;
12023b2fbd54SBarry Smith 
12033a40ed3dSBarry Smith   PetscFunctionBegin;
12043b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1205a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1206b6490206SBarry Smith 
12070513a670SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
12080513a670SBarry Smith 
12093a40ed3dSBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1210a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
12113a40ed3dSBarry Smith   }
12122593348eSBarry Smith 
12132593348eSBarry Smith   *A = 0;
1214f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView);
12152593348eSBarry Smith   PLogObjectCreate(B);
1216b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
121744cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
1218f830108cSBarry Smith   PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps));
12191a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
12201a6a6d98SBarry Smith   if (!flg) {
12217fc0212eSBarry Smith     switch (bs) {
12227fc0212eSBarry Smith     case 1:
1223f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1224f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1225f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1226f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
12277fc0212eSBarry Smith       break;
12284eeb42bcSBarry Smith     case 2:
1229f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1230f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1231f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1232f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
12336d84be18SBarry Smith       break;
1234f327dd97SBarry Smith     case 3:
1235f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1236f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1237f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1238f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
12394eeb42bcSBarry Smith       break;
12406d84be18SBarry Smith     case 4:
1241f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1242f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1243f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1244f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
12456d84be18SBarry Smith       break;
12466d84be18SBarry Smith     case 5:
1247f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1248f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1249f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1250f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
12516d84be18SBarry Smith       break;
125248e9ddb2SSatish Balay     case 7:
1253f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1254f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1255f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
125648e9ddb2SSatish Balay       break;
12577fc0212eSBarry Smith     }
12581a6a6d98SBarry Smith   }
1259b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
1260b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
12612593348eSBarry Smith   B->factor           = 0;
12622593348eSBarry Smith   B->lupivotthreshold = 1.0;
126390f02eecSBarry Smith   B->mapping          = 0;
12642593348eSBarry Smith   b->row              = 0;
12652593348eSBarry Smith   b->col              = 0;
1266e51c0b9cSSatish Balay   b->icol             = 0;
12672593348eSBarry Smith   b->reallocs         = 0;
12682593348eSBarry Smith 
126944cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
127044cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1271b6490206SBarry Smith   b->mbs     = mbs;
1272f2501298SSatish Balay   b->nbs     = nbs;
1273b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
12742593348eSBarry Smith   if (nnz == PETSC_NULL) {
1275119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
12762593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1277b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1278b6490206SBarry Smith     nz = nz*mbs;
12793a40ed3dSBarry Smith   } else {
12802593348eSBarry Smith     nz = 0;
1281b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
12822593348eSBarry Smith   }
12832593348eSBarry Smith 
12842593348eSBarry Smith   /* allocate the matrix space */
12857e67e3f9SSatish Balay   len   = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int);
12862593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
12877e67e3f9SSatish Balay   PetscMemzero(b->a,nz*bs2*sizeof(Scalar));
12887e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
12892593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
12902593348eSBarry Smith   b->i  = b->j + nz;
12912593348eSBarry Smith   b->singlemalloc = 1;
12922593348eSBarry Smith 
1293de6a44a3SBarry Smith   b->i[0] = 0;
1294b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
12952593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
12962593348eSBarry Smith   }
12972593348eSBarry Smith 
1298b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1299b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1300f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1301b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
13022593348eSBarry Smith 
1303b6490206SBarry Smith   b->bs               = bs;
13049df24120SSatish Balay   b->bs2              = bs2;
1305b6490206SBarry Smith   b->mbs              = mbs;
13062593348eSBarry Smith   b->nz               = 0;
130718e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
13082593348eSBarry Smith   b->sorted           = 0;
13092593348eSBarry Smith   b->roworiented      = 1;
13102593348eSBarry Smith   b->nonew            = 0;
13112593348eSBarry Smith   b->diag             = 0;
13122593348eSBarry Smith   b->solve_work       = 0;
1313de6a44a3SBarry Smith   b->mult_work        = 0;
13142593348eSBarry Smith   b->spptr            = 0;
13154e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
13164e220ebcSLois Curfman McInnes 
13172593348eSBarry Smith   *A = B;
13182593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
13192593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
13203a40ed3dSBarry Smith   PetscFunctionReturn(0);
13212593348eSBarry Smith }
13222593348eSBarry Smith 
13235615d1e5SSatish Balay #undef __FUNC__
13245615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ"
1325b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
13262593348eSBarry Smith {
13272593348eSBarry Smith   Mat         C;
1328b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
13299df24120SSatish Balay   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
1330de6a44a3SBarry Smith 
13313a40ed3dSBarry Smith   PetscFunctionBegin;
1332a8c6a408SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
13332593348eSBarry Smith 
13342593348eSBarry Smith   *B = 0;
1335f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView);
13362593348eSBarry Smith   PLogObjectCreate(C);
1337b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
1338c3c66ee5SSatish Balay   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1339b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
1340b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
13412593348eSBarry Smith   C->factor     = A->factor;
13422593348eSBarry Smith   c->row        = 0;
13432593348eSBarry Smith   c->col        = 0;
13442593348eSBarry Smith   C->assembled  = PETSC_TRUE;
13452593348eSBarry Smith 
134644cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
134744cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
134844cd7ae7SLois Curfman McInnes   C->M          = a->m;
134944cd7ae7SLois Curfman McInnes   C->N          = a->n;
135044cd7ae7SLois Curfman McInnes 
1351b6490206SBarry Smith   c->bs         = a->bs;
13529df24120SSatish Balay   c->bs2        = a->bs2;
1353b6490206SBarry Smith   c->mbs        = a->mbs;
1354b6490206SBarry Smith   c->nbs        = a->nbs;
13552593348eSBarry Smith 
1356b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1357b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1358b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
13592593348eSBarry Smith     c->imax[i] = a->imax[i];
13602593348eSBarry Smith     c->ilen[i] = a->ilen[i];
13612593348eSBarry Smith   }
13622593348eSBarry Smith 
13632593348eSBarry Smith   /* allocate the matrix space */
13642593348eSBarry Smith   c->singlemalloc = 1;
13657e67e3f9SSatish Balay   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int));
13662593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
13677e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1368de6a44a3SBarry Smith   c->i  = c->j + nz;
1369b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1370b6490206SBarry Smith   if (mbs > 0) {
1371de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
13722593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
13737e67e3f9SSatish Balay       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar));
13742593348eSBarry Smith     }
13752593348eSBarry Smith   }
13762593348eSBarry Smith 
1377f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
13782593348eSBarry Smith   c->sorted      = a->sorted;
13792593348eSBarry Smith   c->roworiented = a->roworiented;
13802593348eSBarry Smith   c->nonew       = a->nonew;
13812593348eSBarry Smith 
13822593348eSBarry Smith   if (a->diag) {
1383b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1384b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1385b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
13862593348eSBarry Smith       c->diag[i] = a->diag[i];
13872593348eSBarry Smith     }
13882593348eSBarry Smith   }
13892593348eSBarry Smith   else c->diag          = 0;
13902593348eSBarry Smith   c->nz                 = a->nz;
13912593348eSBarry Smith   c->maxnz              = a->maxnz;
13922593348eSBarry Smith   c->solve_work         = 0;
13932593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
13947fc0212eSBarry Smith   c->mult_work          = 0;
13952593348eSBarry Smith   *B = C;
13963a40ed3dSBarry Smith   PetscFunctionReturn(0);
13972593348eSBarry Smith }
13982593348eSBarry Smith 
13995615d1e5SSatish Balay #undef __FUNC__
14005615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ"
140119bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
14022593348eSBarry Smith {
1403b6490206SBarry Smith   Mat_SeqBAIJ  *a;
14042593348eSBarry Smith   Mat          B;
1405de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1406b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
140735aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1408a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1409b6490206SBarry Smith   Scalar       *aa;
141019bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
14112593348eSBarry Smith 
14123a40ed3dSBarry Smith   PetscFunctionBegin;
14131a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1414de6a44a3SBarry Smith   bs2  = bs*bs;
1415b6490206SBarry Smith 
14162593348eSBarry Smith   MPI_Comm_size(comm,&size);
1417cc0fae11SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
141890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
14190752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
1420a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
14212593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
14222593348eSBarry Smith 
1423d64ed03dSBarry Smith   if (header[3] < 0) {
1424a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1425d64ed03dSBarry Smith   }
1426d64ed03dSBarry Smith 
1427a8c6a408SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
142835aab85fSBarry Smith 
142935aab85fSBarry Smith   /*
143035aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
143135aab85fSBarry Smith     divisible by the blocksize
143235aab85fSBarry Smith   */
1433b6490206SBarry Smith   mbs        = M/bs;
143435aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
143535aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
143635aab85fSBarry Smith   else                  mbs++;
143735aab85fSBarry Smith   if (extra_rows) {
1438537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
143935aab85fSBarry Smith   }
1440b6490206SBarry Smith 
14412593348eSBarry Smith   /* read in row lengths */
144235aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
14430752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
144435aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
14452593348eSBarry Smith 
1446b6490206SBarry Smith   /* read in column indices */
144735aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
14480752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr);
144935aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1450b6490206SBarry Smith 
1451b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1452b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1453b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
145435aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
145535aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
145635aab85fSBarry Smith   masked = mask + mbs;
1457b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1458b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
145935aab85fSBarry Smith     nmask = 0;
1460b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1461b6490206SBarry Smith       kmax = rowlengths[rowcount];
1462b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
146335aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
146435aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1465b6490206SBarry Smith       }
1466b6490206SBarry Smith       rowcount++;
1467b6490206SBarry Smith     }
146835aab85fSBarry Smith     browlengths[i] += nmask;
146935aab85fSBarry Smith     /* zero out the mask elements we set */
147035aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1471b6490206SBarry Smith   }
1472b6490206SBarry Smith 
14732593348eSBarry Smith   /* create our matrix */
14743eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
14752593348eSBarry Smith   B = *A;
1476b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
14772593348eSBarry Smith 
14782593348eSBarry Smith   /* set matrix "i" values */
1479de6a44a3SBarry Smith   a->i[0] = 0;
1480b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1481b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1482b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
14832593348eSBarry Smith   }
1484b6490206SBarry Smith   a->nz         = 0;
1485b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
14862593348eSBarry Smith 
1487b6490206SBarry Smith   /* read in nonzero values */
148835aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
14890752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr);
149035aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1491b6490206SBarry Smith 
1492b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1493b6490206SBarry Smith   nzcount = 0; jcount = 0;
1494b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1495b6490206SBarry Smith     nzcountb = nzcount;
149635aab85fSBarry Smith     nmask    = 0;
1497b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1498b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1499b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
150035aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
150135aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1502b6490206SBarry Smith       }
1503b6490206SBarry Smith       rowcount++;
1504b6490206SBarry Smith     }
1505de6a44a3SBarry Smith     /* sort the masked values */
150677c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1507de6a44a3SBarry Smith 
1508b6490206SBarry Smith     /* set "j" values into matrix */
1509b6490206SBarry Smith     maskcount = 1;
151035aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
151135aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1512de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1513b6490206SBarry Smith     }
1514b6490206SBarry Smith     /* set "a" values into matrix */
1515de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1516b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1517b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1518b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1519de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
1520de6a44a3SBarry Smith         block  = mask[tmp] - 1;
1521de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
1522de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
1523b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1524b6490206SBarry Smith       }
1525b6490206SBarry Smith     }
152635aab85fSBarry Smith     /* zero out the mask elements we set */
152735aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1528b6490206SBarry Smith   }
1529a8c6a408SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1530b6490206SBarry Smith 
1531b6490206SBarry Smith   PetscFree(rowlengths);
1532b6490206SBarry Smith   PetscFree(browlengths);
1533b6490206SBarry Smith   PetscFree(aa);
1534b6490206SBarry Smith   PetscFree(jj);
1535b6490206SBarry Smith   PetscFree(mask);
1536b6490206SBarry Smith 
1537b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1538de6a44a3SBarry Smith 
15399c01be13SBarry Smith   ierr = MatView_Private(B); CHKERRQ(ierr);
15403a40ed3dSBarry Smith   PetscFunctionReturn(0);
15412593348eSBarry Smith }
15422593348eSBarry Smith 
15432593348eSBarry Smith 
15442593348eSBarry Smith 
1545