xref: /petsc/src/mat/impls/baij/seq/baij.c (revision b51ba29f2ce68dcf965c343af6f916e8c56c2d8b)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*b51ba29fSSatish Balay static char vcid[] = "$Id: baij.c,v 1.124 1998/01/27 16:11:14 bsmith 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;
1092d61bbb3SSatish Balay 
1102d61bbb3SSatish Balay #if defined(USE_PETSC_LOG)
1112d61bbb3SSatish Balay   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
1122d61bbb3SSatish Balay #endif
1132d61bbb3SSatish Balay   PetscFree(a->a);
1142d61bbb3SSatish Balay   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
1152d61bbb3SSatish Balay   if (a->diag) PetscFree(a->diag);
1162d61bbb3SSatish Balay   if (a->ilen) PetscFree(a->ilen);
1172d61bbb3SSatish Balay   if (a->imax) PetscFree(a->imax);
1182d61bbb3SSatish Balay   if (a->solve_work) PetscFree(a->solve_work);
1192d61bbb3SSatish Balay   if (a->mult_work) PetscFree(a->mult_work);
1202d61bbb3SSatish Balay   PetscFree(a);
1212d61bbb3SSatish Balay   PLogObjectDestroy(A);
1222d61bbb3SSatish Balay   PetscHeaderDestroy(A);
1232d61bbb3SSatish Balay   PetscFunctionReturn(0);
1242d61bbb3SSatish Balay }
1252d61bbb3SSatish Balay 
1262d61bbb3SSatish Balay #undef __FUNC__
1272d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ"
1282d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1292d61bbb3SSatish Balay {
1302d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1312d61bbb3SSatish Balay 
1322d61bbb3SSatish Balay   PetscFunctionBegin;
1332d61bbb3SSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
1342d61bbb3SSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
1352d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
1362d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
1372d61bbb3SSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
1382d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
1392d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
1402d61bbb3SSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
1412d61bbb3SSatish Balay   else if (op == MAT_ROWS_SORTED ||
1422d61bbb3SSatish Balay            op == MAT_ROWS_UNSORTED ||
1432d61bbb3SSatish Balay            op == MAT_SYMMETRIC ||
1442d61bbb3SSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
1452d61bbb3SSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
146*b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
147*b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE) {
148981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
1492d61bbb3SSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
1502d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
1512d61bbb3SSatish Balay   } else {
1522d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
1532d61bbb3SSatish Balay   }
1542d61bbb3SSatish Balay   PetscFunctionReturn(0);
1552d61bbb3SSatish Balay }
1562d61bbb3SSatish Balay 
1572d61bbb3SSatish Balay 
1582d61bbb3SSatish Balay #undef __FUNC__
1592d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ"
1602d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
1612d61bbb3SSatish Balay {
1622d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1632d61bbb3SSatish Balay 
1642d61bbb3SSatish Balay   PetscFunctionBegin;
1652d61bbb3SSatish Balay   if (m) *m = a->m;
1662d61bbb3SSatish Balay   if (n) *n = a->n;
1672d61bbb3SSatish Balay   PetscFunctionReturn(0);
1682d61bbb3SSatish Balay }
1692d61bbb3SSatish Balay 
1702d61bbb3SSatish Balay #undef __FUNC__
1712d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ"
1722d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
1732d61bbb3SSatish Balay {
1742d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1752d61bbb3SSatish Balay 
1762d61bbb3SSatish Balay   PetscFunctionBegin;
1772d61bbb3SSatish Balay   *m = 0; *n = a->m;
1782d61bbb3SSatish Balay   PetscFunctionReturn(0);
1792d61bbb3SSatish Balay }
1802d61bbb3SSatish Balay 
1812d61bbb3SSatish Balay #undef __FUNC__
1822d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ"
1832d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
1842d61bbb3SSatish Balay {
1852d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1862d61bbb3SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
1872d61bbb3SSatish Balay   Scalar      *aa,*v_i,*aa_i;
1882d61bbb3SSatish Balay 
1892d61bbb3SSatish Balay   PetscFunctionBegin;
1902d61bbb3SSatish Balay   bs  = a->bs;
1912d61bbb3SSatish Balay   ai  = a->i;
1922d61bbb3SSatish Balay   aj  = a->j;
1932d61bbb3SSatish Balay   aa  = a->a;
1942d61bbb3SSatish Balay   bs2 = a->bs2;
1952d61bbb3SSatish Balay 
1962d61bbb3SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
1972d61bbb3SSatish Balay 
1982d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
1992d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2002d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2012d61bbb3SSatish Balay   *nz = bs*M;
2022d61bbb3SSatish Balay 
2032d61bbb3SSatish Balay   if (v) {
2042d61bbb3SSatish Balay     *v = 0;
2052d61bbb3SSatish Balay     if (*nz) {
2062d61bbb3SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
2072d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2082d61bbb3SSatish Balay         v_i  = *v + i*bs;
2092d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2102d61bbb3SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
2112d61bbb3SSatish Balay       }
2122d61bbb3SSatish Balay     }
2132d61bbb3SSatish Balay   }
2142d61bbb3SSatish Balay 
2152d61bbb3SSatish Balay   if (idx) {
2162d61bbb3SSatish Balay     *idx = 0;
2172d61bbb3SSatish Balay     if (*nz) {
2182d61bbb3SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
2192d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2202d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2212d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2222d61bbb3SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
2232d61bbb3SSatish Balay       }
2242d61bbb3SSatish Balay     }
2252d61bbb3SSatish Balay   }
2262d61bbb3SSatish Balay   PetscFunctionReturn(0);
2272d61bbb3SSatish Balay }
2282d61bbb3SSatish Balay 
2292d61bbb3SSatish Balay #undef __FUNC__
2302d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ"
2312d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2322d61bbb3SSatish Balay {
2332d61bbb3SSatish Balay   PetscFunctionBegin;
2342d61bbb3SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
2352d61bbb3SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
2362d61bbb3SSatish Balay   PetscFunctionReturn(0);
2372d61bbb3SSatish Balay }
2382d61bbb3SSatish Balay 
2392d61bbb3SSatish Balay #undef __FUNC__
2402d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ"
2412d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
2422d61bbb3SSatish Balay {
2432d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
2442d61bbb3SSatish Balay   Mat         C;
2452d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
2462d61bbb3SSatish Balay   int         *rows,*cols,bs2=a->bs2;
2472d61bbb3SSatish Balay   Scalar      *array=a->a;
2482d61bbb3SSatish Balay 
2492d61bbb3SSatish Balay   PetscFunctionBegin;
2502d61bbb3SSatish Balay   if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
2512d61bbb3SSatish Balay   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
2522d61bbb3SSatish Balay   PetscMemzero(col,(1+nbs)*sizeof(int));
2532d61bbb3SSatish Balay 
2542d61bbb3SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
2552d61bbb3SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
2562d61bbb3SSatish Balay   PetscFree(col);
2572d61bbb3SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
2582d61bbb3SSatish Balay   cols = rows + bs;
2592d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2602d61bbb3SSatish Balay     cols[0] = i*bs;
2612d61bbb3SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
2622d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
2632d61bbb3SSatish Balay     for ( j=0; j<len; j++ ) {
2642d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
2652d61bbb3SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
2662d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
2672d61bbb3SSatish Balay       array += bs2;
2682d61bbb3SSatish Balay     }
2692d61bbb3SSatish Balay   }
2702d61bbb3SSatish Balay   PetscFree(rows);
2712d61bbb3SSatish Balay 
2722d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2732d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2742d61bbb3SSatish Balay 
2752d61bbb3SSatish Balay   if (B != PETSC_NULL) {
2762d61bbb3SSatish Balay     *B = C;
2772d61bbb3SSatish Balay   } else {
2782d61bbb3SSatish Balay     /* This isn't really an in-place transpose */
2792d61bbb3SSatish Balay     PetscFree(a->a);
2802d61bbb3SSatish Balay     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
2812d61bbb3SSatish Balay     if (a->diag) PetscFree(a->diag);
2822d61bbb3SSatish Balay     if (a->ilen) PetscFree(a->ilen);
2832d61bbb3SSatish Balay     if (a->imax) PetscFree(a->imax);
2842d61bbb3SSatish Balay     if (a->solve_work) PetscFree(a->solve_work);
2852d61bbb3SSatish Balay     PetscFree(a);
2862d61bbb3SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
2872d61bbb3SSatish Balay     PetscHeaderDestroy(C);
2882d61bbb3SSatish Balay   }
2892d61bbb3SSatish Balay   PetscFunctionReturn(0);
2902d61bbb3SSatish Balay }
2912d61bbb3SSatish Balay 
2922d61bbb3SSatish Balay 
2932d61bbb3SSatish Balay 
2943b2fbd54SBarry Smith 
2955615d1e5SSatish Balay #undef __FUNC__
296d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary"
297b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
2982593348eSBarry Smith {
299b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3009df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
301b6490206SBarry Smith   Scalar      *aa;
3022593348eSBarry Smith 
3033a40ed3dSBarry Smith   PetscFunctionBegin;
30490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
3052593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
3062593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3073b2fbd54SBarry Smith 
3082593348eSBarry Smith   col_lens[1] = a->m;
3092593348eSBarry Smith   col_lens[2] = a->n;
3107e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3112593348eSBarry Smith 
3122593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
313b6490206SBarry Smith   count = 0;
314b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
315b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
316b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
317b6490206SBarry Smith     }
3182593348eSBarry Smith   }
3190752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3202593348eSBarry Smith   PetscFree(col_lens);
3212593348eSBarry Smith 
3222593348eSBarry Smith   /* store column indices (zero start index) */
32366963ce1SSatish Balay   jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj);
324b6490206SBarry Smith   count = 0;
325b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
326b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
327b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
328b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
329b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3302593348eSBarry Smith         }
3312593348eSBarry Smith       }
332b6490206SBarry Smith     }
333b6490206SBarry Smith   }
3340752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr);
335b6490206SBarry Smith   PetscFree(jj);
3362593348eSBarry Smith 
3372593348eSBarry Smith   /* store nonzero values */
33866963ce1SSatish Balay   aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa);
339b6490206SBarry Smith   count = 0;
340b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
341b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
342b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
343b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
3447e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
345b6490206SBarry Smith         }
346b6490206SBarry Smith       }
347b6490206SBarry Smith     }
348b6490206SBarry Smith   }
3490752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
350b6490206SBarry Smith   PetscFree(aa);
3513a40ed3dSBarry Smith   PetscFunctionReturn(0);
3522593348eSBarry Smith }
3532593348eSBarry Smith 
3545615d1e5SSatish Balay #undef __FUNC__
355d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII"
356b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
3572593348eSBarry Smith {
358b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3599df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
3602593348eSBarry Smith   FILE        *fd;
3612593348eSBarry Smith   char        *outputname;
3622593348eSBarry Smith 
3633a40ed3dSBarry Smith   PetscFunctionBegin;
36490ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
3652593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
36690ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
367639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
3687ddc982cSLois Curfman McInnes     fprintf(fd,"  block size is %d\n",bs);
3692593348eSBarry Smith   }
370639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
371a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
3722593348eSBarry Smith   }
373639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
37444cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
37544cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
37644cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
37744cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
37844cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
3793a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
380766eeae4SLois Curfman McInnes           if (imag(a->a[bs2*k + l*bs + j]) > 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0)
38144cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
38244cd7ae7SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
383766eeae4SLois Curfman McInnes           else if (imag(a->a[bs2*k + l*bs + j]) < 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0)
384766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
385766eeae4SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j]));
38644cd7ae7SLois Curfman McInnes           else if (real(a->a[bs2*k + l*bs + j]) != 0.0)
38744cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
38844cd7ae7SLois Curfman McInnes #else
38944cd7ae7SLois Curfman McInnes           if (a->a[bs2*k + l*bs + j] != 0.0)
39044cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
39144cd7ae7SLois Curfman McInnes #endif
39244cd7ae7SLois Curfman McInnes           }
39344cd7ae7SLois Curfman McInnes         }
39444cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
39544cd7ae7SLois Curfman McInnes       }
39644cd7ae7SLois Curfman McInnes     }
39744cd7ae7SLois Curfman McInnes   }
3982593348eSBarry Smith   else {
399b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
400b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
401b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
402b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
403b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
4043a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
405766eeae4SLois Curfman McInnes           if (imag(a->a[bs2*k + l*bs + j]) > 0.0) {
40688685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
4077e67e3f9SSatish Balay               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
40888685aaeSLois Curfman McInnes           }
409766eeae4SLois Curfman McInnes           else if (imag(a->a[bs2*k + l*bs + j]) < 0.0) {
410766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
411766eeae4SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j]));
412766eeae4SLois Curfman McInnes           }
41388685aaeSLois Curfman McInnes           else {
4147e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
41588685aaeSLois Curfman McInnes           }
41688685aaeSLois Curfman McInnes #else
4177e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
41888685aaeSLois Curfman McInnes #endif
4192593348eSBarry Smith           }
4202593348eSBarry Smith         }
4212593348eSBarry Smith         fprintf(fd,"\n");
4222593348eSBarry Smith       }
4232593348eSBarry Smith     }
424b6490206SBarry Smith   }
4252593348eSBarry Smith   fflush(fd);
4263a40ed3dSBarry Smith   PetscFunctionReturn(0);
4272593348eSBarry Smith }
4282593348eSBarry Smith 
4295615d1e5SSatish Balay #undef __FUNC__
430d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw"
4313270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
4323270192aSSatish Balay {
4333270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
4343270192aSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2;
4353270192aSSatish Balay   double       xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
4363270192aSSatish Balay   Scalar       *aa;
4373270192aSSatish Balay   Draw         draw;
4383270192aSSatish Balay   DrawButton   button;
4393270192aSSatish Balay   PetscTruth   isnull;
4403270192aSSatish Balay 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
4423a40ed3dSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr);
4433a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
4443270192aSSatish Balay 
4453270192aSSatish Balay   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
4463270192aSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
4473270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
4483270192aSSatish Balay   /* loop over matrix elements drawing boxes */
4493270192aSSatish Balay   color = DRAW_BLUE;
4503270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4513270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4523270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4533270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4543270192aSSatish Balay       aa = a->a + j*bs2;
4553270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4563270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4573270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
458b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4593270192aSSatish Balay         }
4603270192aSSatish Balay       }
4613270192aSSatish Balay     }
4623270192aSSatish Balay   }
4633270192aSSatish Balay   color = DRAW_CYAN;
4643270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4653270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4663270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4673270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4683270192aSSatish Balay       aa = a->a + j*bs2;
4693270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4703270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4713270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
472b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4733270192aSSatish Balay         }
4743270192aSSatish Balay       }
4753270192aSSatish Balay     }
4763270192aSSatish Balay   }
4773270192aSSatish Balay 
4783270192aSSatish Balay   color = DRAW_RED;
4793270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4803270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4813270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4823270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4833270192aSSatish Balay       aa = a->a + j*bs2;
4843270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4853270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4863270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
487b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4883270192aSSatish Balay         }
4893270192aSSatish Balay       }
4903270192aSSatish Balay     }
4913270192aSSatish Balay   }
4923270192aSSatish Balay 
49355843e3eSBarry Smith   DrawSynchronizedFlush(draw);
4943270192aSSatish Balay   DrawGetPause(draw,&pause);
4953a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
4963270192aSSatish Balay 
4973270192aSSatish Balay   /* allow the matrix to zoom or shrink */
49855843e3eSBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
4993270192aSSatish Balay   while (button != BUTTON_RIGHT) {
50055843e3eSBarry Smith     DrawSynchronizedClear(draw);
5013270192aSSatish Balay     if (button == BUTTON_LEFT) scale = .5;
5023270192aSSatish Balay     else if (button == BUTTON_CENTER) scale = 2.;
5033270192aSSatish Balay     xl = scale*(xl + w - xc) + xc - w*scale;
5043270192aSSatish Balay     xr = scale*(xr - w - xc) + xc + w*scale;
5053270192aSSatish Balay     yl = scale*(yl + h - yc) + yc - h*scale;
5063270192aSSatish Balay     yr = scale*(yr - h - yc) + yc + h*scale;
5073270192aSSatish Balay     w *= scale; h *= scale;
5083270192aSSatish Balay     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5093270192aSSatish Balay     color = DRAW_BLUE;
5103270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5113270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5123270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5133270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5143270192aSSatish Balay         aa = a->a + j*bs2;
5153270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5163270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5173270192aSSatish Balay             if (PetscReal(*aa++) >=  0.) continue;
518b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5193270192aSSatish Balay           }
5203270192aSSatish Balay         }
5213270192aSSatish Balay       }
5223270192aSSatish Balay     }
5233270192aSSatish Balay     color = DRAW_CYAN;
5243270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5253270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5263270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5273270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5283270192aSSatish Balay         aa = a->a + j*bs2;
5293270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5303270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5313270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
532b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5333270192aSSatish Balay           }
5343270192aSSatish Balay         }
5353270192aSSatish Balay       }
5363270192aSSatish Balay     }
5373270192aSSatish Balay 
5383270192aSSatish Balay     color = DRAW_RED;
5393270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5403270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5413270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5423270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5433270192aSSatish Balay         aa = a->a + j*bs2;
5443270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5453270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5463270192aSSatish Balay             if (PetscReal(*aa++) <= 0.) continue;
547b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5483270192aSSatish Balay           }
5493270192aSSatish Balay         }
5503270192aSSatish Balay       }
5513270192aSSatish Balay     }
55255843e3eSBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
5533270192aSSatish Balay   }
5543a40ed3dSBarry Smith   PetscFunctionReturn(0);
5553270192aSSatish Balay }
5563270192aSSatish Balay 
5575615d1e5SSatish Balay #undef __FUNC__
558d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ"
5598f6be9afSLois Curfman McInnes int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
5602593348eSBarry Smith {
5612593348eSBarry Smith   Mat         A = (Mat) obj;
56219bcc07fSBarry Smith   ViewerType  vtype;
56319bcc07fSBarry Smith   int         ierr;
5642593348eSBarry Smith 
5653a40ed3dSBarry Smith   PetscFunctionBegin;
56619bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
56719bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
568a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported");
5693a40ed3dSBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
5703a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
5713a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
5723a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
5733a40ed3dSBarry Smith   } else if (vtype == DRAW_VIEWER) {
5743a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
5752593348eSBarry Smith   }
5763a40ed3dSBarry Smith   PetscFunctionReturn(0);
5772593348eSBarry Smith }
578b6490206SBarry Smith 
579cd0e1443SSatish Balay 
5805615d1e5SSatish Balay #undef __FUNC__
5812d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ"
5822d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
583cd0e1443SSatish Balay {
584cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5852d61bbb3SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
5862d61bbb3SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
5872d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
5882d61bbb3SSatish Balay   Scalar     *ap, *aa = a->a, zero = 0.0;
589cd0e1443SSatish Balay 
5903a40ed3dSBarry Smith   PetscFunctionBegin;
5912d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
592cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
593a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
594a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
5952d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
5962c3acbe9SBarry Smith     nrow = ailen[brow];
5972d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
598a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
599a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6002d61bbb3SSatish Balay       col  = in[l] ;
6012d61bbb3SSatish Balay       bcol = col/bs;
6022d61bbb3SSatish Balay       cidx = col%bs;
6032d61bbb3SSatish Balay       ridx = row%bs;
6042d61bbb3SSatish Balay       high = nrow;
6052d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6062d61bbb3SSatish Balay       while (high-low > 5) {
607cd0e1443SSatish Balay         t = (low+high)/2;
608cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
609cd0e1443SSatish Balay         else             low  = t;
610cd0e1443SSatish Balay       }
611cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
612cd0e1443SSatish Balay         if (rp[i] > bcol) break;
613cd0e1443SSatish Balay         if (rp[i] == bcol) {
6142d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6152d61bbb3SSatish Balay           goto finished;
616cd0e1443SSatish Balay         }
617cd0e1443SSatish Balay       }
6182d61bbb3SSatish Balay       *v++ = zero;
6192d61bbb3SSatish Balay       finished:;
620cd0e1443SSatish Balay     }
621cd0e1443SSatish Balay   }
6223a40ed3dSBarry Smith   PetscFunctionReturn(0);
623cd0e1443SSatish Balay }
624cd0e1443SSatish Balay 
6252d61bbb3SSatish Balay 
6265615d1e5SSatish Balay #undef __FUNC__
62705a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ"
62892c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
62992c4ed94SBarry Smith {
63092c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6318a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
632dd9472c6SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented;
633dd9472c6SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
6340e324ae4SSatish Balay   Scalar      *ap,*value=v,*aa=a->a,*bap;
63592c4ed94SBarry Smith 
6363a40ed3dSBarry Smith   PetscFunctionBegin;
6370e324ae4SSatish Balay   if (roworiented) {
6380e324ae4SSatish Balay     stepval = (n-1)*bs;
6390e324ae4SSatish Balay   } else {
6400e324ae4SSatish Balay     stepval = (m-1)*bs;
6410e324ae4SSatish Balay   }
64292c4ed94SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
64392c4ed94SBarry Smith     row  = im[k];
6443a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
645a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
646a8c6a408SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
64792c4ed94SBarry Smith #endif
64892c4ed94SBarry Smith     rp   = aj + ai[row];
64992c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
65092c4ed94SBarry Smith     rmax = imax[row];
65192c4ed94SBarry Smith     nrow = ailen[row];
65292c4ed94SBarry Smith     low  = 0;
65392c4ed94SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
6543a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
655a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
656a8c6a408SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
65792c4ed94SBarry Smith #endif
65892c4ed94SBarry Smith       col = in[l];
65992c4ed94SBarry Smith       if (roworiented) {
6600e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
6610e324ae4SSatish Balay       } else {
6620e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
66392c4ed94SBarry Smith       }
66492c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
66592c4ed94SBarry Smith       while (high-low > 7) {
66692c4ed94SBarry Smith         t = (low+high)/2;
66792c4ed94SBarry Smith         if (rp[t] > col) high = t;
66892c4ed94SBarry Smith         else             low  = t;
66992c4ed94SBarry Smith       }
67092c4ed94SBarry Smith       for ( i=low; i<high; i++ ) {
67192c4ed94SBarry Smith         if (rp[i] > col) break;
67292c4ed94SBarry Smith         if (rp[i] == col) {
6738a84c255SSatish Balay           bap  = ap +  bs2*i;
6740e324ae4SSatish Balay           if (roworiented) {
6758a84c255SSatish Balay             if (is == ADD_VALUES) {
676dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
677dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
6788a84c255SSatish Balay                   bap[jj] += *value++;
679dd9472c6SBarry Smith                 }
680dd9472c6SBarry Smith               }
6810e324ae4SSatish Balay             } else {
682dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
683dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
6840e324ae4SSatish Balay                   bap[jj] = *value++;
6858a84c255SSatish Balay                 }
686dd9472c6SBarry Smith               }
687dd9472c6SBarry Smith             }
6880e324ae4SSatish Balay           } else {
6890e324ae4SSatish Balay             if (is == ADD_VALUES) {
690dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
691dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
6920e324ae4SSatish Balay                   *bap++ += *value++;
693dd9472c6SBarry Smith                 }
694dd9472c6SBarry Smith               }
6950e324ae4SSatish Balay             } else {
696dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
697dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
6980e324ae4SSatish Balay                   *bap++  = *value++;
6990e324ae4SSatish Balay                 }
7008a84c255SSatish Balay               }
701dd9472c6SBarry Smith             }
702dd9472c6SBarry Smith           }
703f1241b54SBarry Smith           goto noinsert2;
70492c4ed94SBarry Smith         }
70592c4ed94SBarry Smith       }
70689280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
707a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
70892c4ed94SBarry Smith       if (nrow >= rmax) {
70992c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
71092c4ed94SBarry Smith         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
71192c4ed94SBarry Smith         Scalar *new_a;
71292c4ed94SBarry Smith 
713a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
71489280ab3SLois Curfman McInnes 
71592c4ed94SBarry Smith         /* malloc new storage space */
71692c4ed94SBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
71792c4ed94SBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
71892c4ed94SBarry Smith         new_j   = (int *) (new_a + bs2*new_nz);
71992c4ed94SBarry Smith         new_i   = new_j + new_nz;
72092c4ed94SBarry Smith 
72192c4ed94SBarry Smith         /* copy over old data into new slots */
72292c4ed94SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
72392c4ed94SBarry Smith         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
72492c4ed94SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));
72592c4ed94SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow);
726dd9472c6SBarry Smith         PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));
72792c4ed94SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(Scalar));
72892c4ed94SBarry Smith         PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
729dd9472c6SBarry Smith         PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(Scalar));
73092c4ed94SBarry Smith         /* free up old matrix storage */
73192c4ed94SBarry Smith         PetscFree(a->a);
73292c4ed94SBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
73392c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
73492c4ed94SBarry Smith         a->singlemalloc = 1;
73592c4ed94SBarry Smith 
73692c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
73792c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
73892c4ed94SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
73992c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
74092c4ed94SBarry Smith         a->reallocs++;
74192c4ed94SBarry Smith         a->nz++;
74292c4ed94SBarry Smith       }
74392c4ed94SBarry Smith       N = nrow++ - 1;
74492c4ed94SBarry Smith       /* shift up all the later entries in this row */
74592c4ed94SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
74692c4ed94SBarry Smith         rp[ii+1] = rp[ii];
74792c4ed94SBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
74892c4ed94SBarry Smith       }
74992c4ed94SBarry Smith       if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
75092c4ed94SBarry Smith       rp[i] = col;
7518a84c255SSatish Balay       bap   = ap +  bs2*i;
7520e324ae4SSatish Balay       if (roworiented) {
753dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval) {
754dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs ) {
7550e324ae4SSatish Balay             bap[jj] = *value++;
756dd9472c6SBarry Smith           }
757dd9472c6SBarry Smith         }
7580e324ae4SSatish Balay       } else {
759dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval ) {
760dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++ ) {
7610e324ae4SSatish Balay             *bap++  = *value++;
7620e324ae4SSatish Balay           }
763dd9472c6SBarry Smith         }
764dd9472c6SBarry Smith       }
765f1241b54SBarry Smith       noinsert2:;
76692c4ed94SBarry Smith       low = i;
76792c4ed94SBarry Smith     }
76892c4ed94SBarry Smith     ailen[row] = nrow;
76992c4ed94SBarry Smith   }
7703a40ed3dSBarry Smith   PetscFunctionReturn(0);
77192c4ed94SBarry Smith }
77292c4ed94SBarry Smith 
773f2501298SSatish Balay 
7745615d1e5SSatish Balay #undef __FUNC__
7755615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ"
7768f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
777584200bdSSatish Balay {
778584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
779584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
780a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
78143ee02c3SBarry Smith   int        mbs = a->mbs, bs2 = a->bs2,rmax = 0;
782584200bdSSatish Balay   Scalar     *aa = a->a, *ap;
783584200bdSSatish Balay 
7843a40ed3dSBarry Smith   PetscFunctionBegin;
7853a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
786584200bdSSatish Balay 
78743ee02c3SBarry Smith   if (m) rmax = ailen[0];
788584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
789584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
790584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
791d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
792584200bdSSatish Balay     if (fshift) {
793a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
794584200bdSSatish Balay       N = ailen[i];
795584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
796584200bdSSatish Balay         ip[j-fshift] = ip[j];
7977e67e3f9SSatish Balay         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar));
798584200bdSSatish Balay       }
799584200bdSSatish Balay     }
800584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
801584200bdSSatish Balay   }
802584200bdSSatish Balay   if (mbs) {
803584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
804584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
805584200bdSSatish Balay   }
806584200bdSSatish Balay   /* reset ilen and imax for each row */
807584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
808584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
809584200bdSSatish Balay   }
810a7c10996SSatish Balay   a->nz = ai[mbs];
811584200bdSSatish Balay 
812584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
813584200bdSSatish Balay   if (fshift && a->diag) {
814584200bdSSatish Balay     PetscFree(a->diag);
815584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
816584200bdSSatish Balay     a->diag = 0;
817584200bdSSatish Balay   }
8183d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
819219d9a1aSLois Curfman McInnes            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
8203d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
821584200bdSSatish Balay            a->reallocs);
822d402145bSBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
823e2f3b5e9SSatish Balay   a->reallocs          = 0;
8244e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
8254e220ebcSLois Curfman McInnes 
8263a40ed3dSBarry Smith   PetscFunctionReturn(0);
827584200bdSSatish Balay }
828584200bdSSatish Balay 
8295a838052SSatish Balay 
830d9b7c43dSSatish Balay /* idx should be of length atlease bs */
8315615d1e5SSatish Balay #undef __FUNC__
8325615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block"
833d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg)
834d9b7c43dSSatish Balay {
835d9b7c43dSSatish Balay   int i,row;
8363a40ed3dSBarry Smith 
8373a40ed3dSBarry Smith   PetscFunctionBegin;
838d9b7c43dSSatish Balay   row = idx[0];
8393a40ed3dSBarry Smith   if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
840d9b7c43dSSatish Balay 
841d9b7c43dSSatish Balay   for ( i=1; i<bs; i++ ) {
8423a40ed3dSBarry Smith     if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
843d9b7c43dSSatish Balay   }
844d9b7c43dSSatish Balay   *flg = PETSC_TRUE;
8453a40ed3dSBarry Smith   PetscFunctionReturn(0);
846d9b7c43dSSatish Balay }
847d9b7c43dSSatish Balay 
8485615d1e5SSatish Balay #undef __FUNC__
8495615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ"
8508f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
851d9b7c43dSSatish Balay {
852d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
853d9b7c43dSSatish Balay   IS          is_local;
854d9b7c43dSSatish Balay   int         ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2;
855d9b7c43dSSatish Balay   PetscTruth  flg;
856d9b7c43dSSatish Balay   Scalar      zero = 0.0,*aa;
857d9b7c43dSSatish Balay 
8583a40ed3dSBarry Smith   PetscFunctionBegin;
859d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
860d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
861d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
862537820f0SBarry Smith   ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr);
863d9b7c43dSSatish Balay   ierr = ISSort(is_local); CHKERRQ(ierr);
864d9b7c43dSSatish Balay   ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr);
865d9b7c43dSSatish Balay 
866d9b7c43dSSatish Balay   i = 0;
867d9b7c43dSSatish Balay   while (i < is_n) {
868a8c6a408SBarry Smith     if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
869d9b7c43dSSatish Balay     flg = PETSC_FALSE;
870d9b7c43dSSatish Balay     if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); }
871d9b7c43dSSatish Balay     count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs;
872d9b7c43dSSatish Balay     aa    = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs);
8732d61bbb3SSatish Balay     if (flg) { /* There exists a block of rows to be Zerowed */
874a07cd24cSSatish Balay       if (baij->ilen[rows[i]/bs] > 0) {
8752d61bbb3SSatish Balay         PetscMemzero(aa,count*bs*sizeof(Scalar));
8762d61bbb3SSatish Balay         baij->ilen[rows[i]/bs] = 1;
8772d61bbb3SSatish Balay         baij->j[baij->i[rows[i]/bs]] = rows[i]/bs;
878a07cd24cSSatish Balay       }
8792d61bbb3SSatish Balay       i += bs;
8802d61bbb3SSatish Balay     } else { /* Zero out only the requested row */
881d9b7c43dSSatish Balay       for ( j=0; j<count; j++ ) {
882d9b7c43dSSatish Balay         aa[0] = zero;
883d9b7c43dSSatish Balay         aa+=bs;
884d9b7c43dSSatish Balay       }
885d9b7c43dSSatish Balay       i++;
886d9b7c43dSSatish Balay     }
887d9b7c43dSSatish Balay   }
888d9b7c43dSSatish Balay   if (diag) {
889d9b7c43dSSatish Balay     for ( j=0; j<is_n; j++ ) {
890d9b7c43dSSatish Balay       ierr = (*A->ops.setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
8912d61bbb3SSatish Balay       /* ierr = MatSetValues(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); */
892d9b7c43dSSatish Balay     }
893d9b7c43dSSatish Balay   }
894d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr);
895d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr);
896d9b7c43dSSatish Balay   ierr = ISDestroy(is_local); CHKERRQ(ierr);
8979a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
8983a40ed3dSBarry Smith   PetscFunctionReturn(0);
899d9b7c43dSSatish Balay }
9001c351548SSatish Balay 
9015615d1e5SSatish Balay #undef __FUNC__
9022d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ"
9032d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
9042d61bbb3SSatish Balay {
9052d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9062d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
9072d61bbb3SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
9082d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
9092d61bbb3SSatish Balay   int          ridx,cidx,bs2=a->bs2;
9102d61bbb3SSatish Balay   Scalar      *ap,value,*aa=a->a,*bap;
9112d61bbb3SSatish Balay 
9122d61bbb3SSatish Balay   PetscFunctionBegin;
9132d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
9142d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
9152d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9162d61bbb3SSatish Balay     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
9172d61bbb3SSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
9182d61bbb3SSatish Balay #endif
9192d61bbb3SSatish Balay     rp   = aj + ai[brow];
9202d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
9212d61bbb3SSatish Balay     rmax = imax[brow];
9222d61bbb3SSatish Balay     nrow = ailen[brow];
9232d61bbb3SSatish Balay     low  = 0;
9242d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
9252d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9262d61bbb3SSatish Balay       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
9272d61bbb3SSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
9282d61bbb3SSatish Balay #endif
9292d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
9302d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
9312d61bbb3SSatish Balay       if (roworiented) {
9322d61bbb3SSatish Balay         value = *v++;
9332d61bbb3SSatish Balay       } else {
9342d61bbb3SSatish Balay         value = v[k + l*m];
9352d61bbb3SSatish Balay       }
9362d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
9372d61bbb3SSatish Balay       while (high-low > 7) {
9382d61bbb3SSatish Balay         t = (low+high)/2;
9392d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
9402d61bbb3SSatish Balay         else              low  = t;
9412d61bbb3SSatish Balay       }
9422d61bbb3SSatish Balay       for ( i=low; i<high; i++ ) {
9432d61bbb3SSatish Balay         if (rp[i] > bcol) break;
9442d61bbb3SSatish Balay         if (rp[i] == bcol) {
9452d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
9462d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
9472d61bbb3SSatish Balay           else                  *bap  = value;
9482d61bbb3SSatish Balay           goto noinsert1;
9492d61bbb3SSatish Balay         }
9502d61bbb3SSatish Balay       }
9512d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
9522d61bbb3SSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9532d61bbb3SSatish Balay       if (nrow >= rmax) {
9542d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
9552d61bbb3SSatish Balay         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
9562d61bbb3SSatish Balay         Scalar *new_a;
9572d61bbb3SSatish Balay 
9582d61bbb3SSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9592d61bbb3SSatish Balay 
9602d61bbb3SSatish Balay         /* Malloc new storage space */
9612d61bbb3SSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
9622d61bbb3SSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
9632d61bbb3SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
9642d61bbb3SSatish Balay         new_i   = new_j + new_nz;
9652d61bbb3SSatish Balay 
9662d61bbb3SSatish Balay         /* copy over old data into new slots */
9672d61bbb3SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
9682d61bbb3SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
9692d61bbb3SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
9702d61bbb3SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
9712d61bbb3SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
9722d61bbb3SSatish Balay                                                            len*sizeof(int));
9732d61bbb3SSatish Balay         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar));
9742d61bbb3SSatish Balay         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
9752d61bbb3SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
9762d61bbb3SSatish Balay                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar));
9772d61bbb3SSatish Balay         /* free up old matrix storage */
9782d61bbb3SSatish Balay         PetscFree(a->a);
9792d61bbb3SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
9802d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
9812d61bbb3SSatish Balay         a->singlemalloc = 1;
9822d61bbb3SSatish Balay 
9832d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
9842d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
9852d61bbb3SSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
9862d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
9872d61bbb3SSatish Balay         a->reallocs++;
9882d61bbb3SSatish Balay         a->nz++;
9892d61bbb3SSatish Balay       }
9902d61bbb3SSatish Balay       N = nrow++ - 1;
9912d61bbb3SSatish Balay       /* shift up all the later entries in this row */
9922d61bbb3SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
9932d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
9942d61bbb3SSatish Balay         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
9952d61bbb3SSatish Balay       }
9962d61bbb3SSatish Balay       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
9972d61bbb3SSatish Balay       rp[i]                      = bcol;
9982d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
9992d61bbb3SSatish Balay       noinsert1:;
10002d61bbb3SSatish Balay       low = i;
10012d61bbb3SSatish Balay     }
10022d61bbb3SSatish Balay     ailen[brow] = nrow;
10032d61bbb3SSatish Balay   }
10042d61bbb3SSatish Balay   PetscFunctionReturn(0);
10052d61bbb3SSatish Balay }
10062d61bbb3SSatish Balay 
10072d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
10082d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
10092d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
1010d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*);
10112d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**);
10122d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec);
10132d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
10142d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat);
10152d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *);
10162d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*);
10172d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec);
10182d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
10192d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
10202d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat);
10212d61bbb3SSatish Balay 
10222d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
10232d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
10242d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
10252d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
10262d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
10272d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
10282d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
10292d61bbb3SSatish Balay 
10302d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
10312d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
10322d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
10332d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
10342d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
10352d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
10362d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
10372d61bbb3SSatish Balay 
10382d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
10392d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
10402d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
10412d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
10422d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
10432d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
10442d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
10452d61bbb3SSatish Balay 
10462d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
10472d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
10482d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
10492d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
10502d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
10512d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
10522d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
10532d61bbb3SSatish Balay 
10542d61bbb3SSatish Balay /*
10552d61bbb3SSatish Balay      note: This can only work for identity for row and col. It would
10562d61bbb3SSatish Balay    be good to check this and otherwise generate an error.
10572d61bbb3SSatish Balay */
10582d61bbb3SSatish Balay #undef __FUNC__
10592d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ"
10602d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
10612d61bbb3SSatish Balay {
10622d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
10632d61bbb3SSatish Balay   Mat         outA;
10642d61bbb3SSatish Balay   int         ierr;
10652d61bbb3SSatish Balay 
10662d61bbb3SSatish Balay   PetscFunctionBegin;
10672d61bbb3SSatish Balay   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
10682d61bbb3SSatish Balay 
10692d61bbb3SSatish Balay   outA          = inA;
10702d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
10712d61bbb3SSatish Balay   a->row        = row;
10722d61bbb3SSatish Balay   a->col        = col;
10732d61bbb3SSatish Balay 
10742d61bbb3SSatish Balay   if (!a->solve_work) {
10752d61bbb3SSatish Balay     a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
10762d61bbb3SSatish Balay     PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
10772d61bbb3SSatish Balay   }
10782d61bbb3SSatish Balay 
10792d61bbb3SSatish Balay   if (!a->diag) {
10802d61bbb3SSatish Balay     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
10812d61bbb3SSatish Balay   }
10822d61bbb3SSatish Balay   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
10832d61bbb3SSatish Balay 
10842d61bbb3SSatish Balay   /*
10852d61bbb3SSatish Balay       Blocksize 4 has a special faster solver for ILU(0) factorization
10862d61bbb3SSatish Balay     with natural ordering
10872d61bbb3SSatish Balay   */
10882d61bbb3SSatish Balay   if (a->bs == 4) {
10892d61bbb3SSatish Balay     inA->ops.solve = MatSolve_SeqBAIJ_4_NaturalOrdering;
10902d61bbb3SSatish Balay   }
10912d61bbb3SSatish Balay 
10922d61bbb3SSatish Balay   PetscFunctionReturn(0);
10932d61bbb3SSatish Balay }
10942d61bbb3SSatish Balay #undef __FUNC__
1095d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ"
1096ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1097ba4ca20aSSatish Balay {
1098ba4ca20aSSatish Balay   static int called = 0;
1099ba4ca20aSSatish Balay   MPI_Comm   comm = A->comm;
1100ba4ca20aSSatish Balay 
11013a40ed3dSBarry Smith   PetscFunctionBegin;
11023a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
110376be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");
110476be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");
11053a40ed3dSBarry Smith   PetscFunctionReturn(0);
1106ba4ca20aSSatish Balay }
1107d9b7c43dSSatish Balay 
11082593348eSBarry Smith /* -------------------------------------------------------------------*/
1109cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ,
11109867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
1111f44d3a6dSSatish Balay        MatMult_SeqBAIJ_N,MatMultAdd_SeqBAIJ_N,
1112faf6580fSSatish Balay        MatMultTrans_SeqBAIJ,MatMultTransAdd_SeqBAIJ,
11131a6a6d98SBarry Smith        MatSolve_SeqBAIJ_N,0,
1114b6490206SBarry Smith        0,0,
1115de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
1116b6490206SBarry Smith        0,
1117f2501298SSatish Balay        MatTranspose_SeqBAIJ,
111817e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
11199867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
1120584200bdSSatish Balay        0,MatAssemblyEnd_SeqBAIJ,
1121b6490206SBarry Smith        0,
1122d9b7c43dSSatish Balay        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,MatZeroRows_SeqBAIJ,
11237fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
1124b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
1125de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
1126d402145bSBarry Smith        0,0,
1127b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
1128b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
1129af015674SSatish Balay        MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ,
11307e67e3f9SSatish Balay        MatGetValues_SeqBAIJ,0,
1131ba4ca20aSSatish Balay        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
11323b2fbd54SBarry Smith        0,0,0,MatGetBlockSize_SeqBAIJ,
11333b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
113492c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
113592c4ed94SBarry Smith        0,0,0,0,0,0,
1136d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1137d3825aa8SBarry Smith        MatGetSubMatrix_SeqBAIJ};
11382593348eSBarry Smith 
11395615d1e5SSatish Balay #undef __FUNC__
11405615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ"
11412593348eSBarry Smith /*@C
114244cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
114344cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
114444cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
11452bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
11462bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
11472593348eSBarry Smith 
11482593348eSBarry Smith    Input Parameters:
1149029af93fSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
1150b6490206SBarry Smith .  bs - size of block
11512593348eSBarry Smith .  m - number of rows
11522593348eSBarry Smith .  n - number of columns
1153b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
11542bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of block nonzeros in the various block rows
11552bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
11562593348eSBarry Smith 
11572593348eSBarry Smith    Output Parameter:
11582593348eSBarry Smith .  A - the matrix
11592593348eSBarry Smith 
11600513a670SBarry Smith    Options Database Keys:
11610513a670SBarry Smith $    -mat_no_unroll - uses code that does not unroll the loops in the
11620effe34fSLois Curfman McInnes $                     block calculations (much slower)
11630513a670SBarry Smith $    -mat_block_size - size of the blocks to use
11640513a670SBarry Smith 
11652593348eSBarry Smith    Notes:
116644cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
11672593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
116844cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
11692593348eSBarry Smith 
11702593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
11712593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
11722593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
11736da5968aSLois Curfman McInnes    matrices.
11742593348eSBarry Smith 
117544cd7ae7SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
11762593348eSBarry Smith @*/
1177b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
11782593348eSBarry Smith {
11792593348eSBarry Smith   Mat         B;
1180b6490206SBarry Smith   Mat_SeqBAIJ *b;
11813b2fbd54SBarry Smith   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size;
11823b2fbd54SBarry Smith 
11833a40ed3dSBarry Smith   PetscFunctionBegin;
11843b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1185a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1186b6490206SBarry Smith 
11870513a670SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
11880513a670SBarry Smith 
11893a40ed3dSBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1190a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
11913a40ed3dSBarry Smith   }
11922593348eSBarry Smith 
11932593348eSBarry Smith   *A = 0;
1194d4bb536fSBarry Smith   PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView);
11952593348eSBarry Smith   PLogObjectCreate(B);
1196b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
119744cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
11982593348eSBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
11991a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
12001a6a6d98SBarry Smith   if (!flg) {
12017fc0212eSBarry Smith     switch (bs) {
12027fc0212eSBarry Smith     case 1:
12037fc0212eSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
12041a6a6d98SBarry Smith       B->ops.solve           = MatSolve_SeqBAIJ_1;
120539b95ed1SBarry Smith       B->ops.mult            = MatMult_SeqBAIJ_1;
1206f44d3a6dSSatish Balay       B->ops.multadd         = MatMultAdd_SeqBAIJ_1;
12077fc0212eSBarry Smith       break;
12084eeb42bcSBarry Smith     case 2:
12094eeb42bcSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
12101a6a6d98SBarry Smith       B->ops.solve           = MatSolve_SeqBAIJ_2;
121139b95ed1SBarry Smith       B->ops.mult            = MatMult_SeqBAIJ_2;
1212f44d3a6dSSatish Balay       B->ops.multadd         = MatMultAdd_SeqBAIJ_2;
12136d84be18SBarry Smith       break;
1214f327dd97SBarry Smith     case 3:
1215f327dd97SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
12161a6a6d98SBarry Smith       B->ops.solve           = MatSolve_SeqBAIJ_3;
121739b95ed1SBarry Smith       B->ops.mult            = MatMult_SeqBAIJ_3;
1218f44d3a6dSSatish Balay       B->ops.multadd         = MatMultAdd_SeqBAIJ_3;
12194eeb42bcSBarry Smith       break;
12206d84be18SBarry Smith     case 4:
12216d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
12221a6a6d98SBarry Smith       B->ops.solve           = MatSolve_SeqBAIJ_4;
122339b95ed1SBarry Smith       B->ops.mult            = MatMult_SeqBAIJ_4;
1224f44d3a6dSSatish Balay       B->ops.multadd         = MatMultAdd_SeqBAIJ_4;
12256d84be18SBarry Smith       break;
12266d84be18SBarry Smith     case 5:
12276d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
12281a6a6d98SBarry Smith       B->ops.solve           = MatSolve_SeqBAIJ_5;
122939b95ed1SBarry Smith       B->ops.mult            = MatMult_SeqBAIJ_5;
1230f44d3a6dSSatish Balay       B->ops.multadd         = MatMultAdd_SeqBAIJ_5;
12316d84be18SBarry Smith       break;
123248e9ddb2SSatish Balay     case 7:
123348e9ddb2SSatish Balay       B->ops.mult            = MatMult_SeqBAIJ_7;
123448e9ddb2SSatish Balay       B->ops.solve           = MatSolve_SeqBAIJ_7;
123548e9ddb2SSatish Balay       B->ops.multadd         = MatMultAdd_SeqBAIJ_7;
123648e9ddb2SSatish Balay       break;
12377fc0212eSBarry Smith     }
12381a6a6d98SBarry Smith   }
1239b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
1240b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
12412593348eSBarry Smith   B->factor           = 0;
12422593348eSBarry Smith   B->lupivotthreshold = 1.0;
124390f02eecSBarry Smith   B->mapping          = 0;
12442593348eSBarry Smith   b->row              = 0;
12452593348eSBarry Smith   b->col              = 0;
12462593348eSBarry Smith   b->reallocs         = 0;
12472593348eSBarry Smith 
124844cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
124944cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1250b6490206SBarry Smith   b->mbs     = mbs;
1251f2501298SSatish Balay   b->nbs     = nbs;
1252b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
12532593348eSBarry Smith   if (nnz == PETSC_NULL) {
1254119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
12552593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1256b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1257b6490206SBarry Smith     nz = nz*mbs;
12583a40ed3dSBarry Smith   } else {
12592593348eSBarry Smith     nz = 0;
1260b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
12612593348eSBarry Smith   }
12622593348eSBarry Smith 
12632593348eSBarry Smith   /* allocate the matrix space */
12647e67e3f9SSatish Balay   len   = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int);
12652593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
12667e67e3f9SSatish Balay   PetscMemzero(b->a,nz*bs2*sizeof(Scalar));
12677e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
12682593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
12692593348eSBarry Smith   b->i  = b->j + nz;
12702593348eSBarry Smith   b->singlemalloc = 1;
12712593348eSBarry Smith 
1272de6a44a3SBarry Smith   b->i[0] = 0;
1273b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
12742593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
12752593348eSBarry Smith   }
12762593348eSBarry Smith 
1277b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1278b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1279f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1280b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
12812593348eSBarry Smith 
1282b6490206SBarry Smith   b->bs               = bs;
12839df24120SSatish Balay   b->bs2              = bs2;
1284b6490206SBarry Smith   b->mbs              = mbs;
12852593348eSBarry Smith   b->nz               = 0;
128618e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
12872593348eSBarry Smith   b->sorted           = 0;
12882593348eSBarry Smith   b->roworiented      = 1;
12892593348eSBarry Smith   b->nonew            = 0;
12902593348eSBarry Smith   b->diag             = 0;
12912593348eSBarry Smith   b->solve_work       = 0;
1292de6a44a3SBarry Smith   b->mult_work        = 0;
12932593348eSBarry Smith   b->spptr            = 0;
12944e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
12954e220ebcSLois Curfman McInnes 
12962593348eSBarry Smith   *A = B;
12972593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
12982593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
12993a40ed3dSBarry Smith   PetscFunctionReturn(0);
13002593348eSBarry Smith }
13012593348eSBarry Smith 
13025615d1e5SSatish Balay #undef __FUNC__
13035615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ"
1304b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
13052593348eSBarry Smith {
13062593348eSBarry Smith   Mat         C;
1307b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
13089df24120SSatish Balay   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
1309de6a44a3SBarry Smith 
13103a40ed3dSBarry Smith   PetscFunctionBegin;
1311a8c6a408SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
13122593348eSBarry Smith 
13132593348eSBarry Smith   *B = 0;
1314d4bb536fSBarry Smith   PetscHeaderCreate(C,_p_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView);
13152593348eSBarry Smith   PLogObjectCreate(C);
1316b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
13172593348eSBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1318b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
1319b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
13202593348eSBarry Smith   C->factor     = A->factor;
13212593348eSBarry Smith   c->row        = 0;
13222593348eSBarry Smith   c->col        = 0;
13232593348eSBarry Smith   C->assembled  = PETSC_TRUE;
13242593348eSBarry Smith 
132544cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
132644cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
132744cd7ae7SLois Curfman McInnes   C->M          = a->m;
132844cd7ae7SLois Curfman McInnes   C->N          = a->n;
132944cd7ae7SLois Curfman McInnes 
1330b6490206SBarry Smith   c->bs         = a->bs;
13319df24120SSatish Balay   c->bs2        = a->bs2;
1332b6490206SBarry Smith   c->mbs        = a->mbs;
1333b6490206SBarry Smith   c->nbs        = a->nbs;
13342593348eSBarry Smith 
1335b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1336b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1337b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
13382593348eSBarry Smith     c->imax[i] = a->imax[i];
13392593348eSBarry Smith     c->ilen[i] = a->ilen[i];
13402593348eSBarry Smith   }
13412593348eSBarry Smith 
13422593348eSBarry Smith   /* allocate the matrix space */
13432593348eSBarry Smith   c->singlemalloc = 1;
13447e67e3f9SSatish Balay   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int));
13452593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
13467e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1347de6a44a3SBarry Smith   c->i  = c->j + nz;
1348b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1349b6490206SBarry Smith   if (mbs > 0) {
1350de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
13512593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
13527e67e3f9SSatish Balay       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar));
13532593348eSBarry Smith     }
13542593348eSBarry Smith   }
13552593348eSBarry Smith 
1356f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
13572593348eSBarry Smith   c->sorted      = a->sorted;
13582593348eSBarry Smith   c->roworiented = a->roworiented;
13592593348eSBarry Smith   c->nonew       = a->nonew;
13602593348eSBarry Smith 
13612593348eSBarry Smith   if (a->diag) {
1362b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1363b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1364b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
13652593348eSBarry Smith       c->diag[i] = a->diag[i];
13662593348eSBarry Smith     }
13672593348eSBarry Smith   }
13682593348eSBarry Smith   else c->diag          = 0;
13692593348eSBarry Smith   c->nz                 = a->nz;
13702593348eSBarry Smith   c->maxnz              = a->maxnz;
13712593348eSBarry Smith   c->solve_work         = 0;
13722593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
13737fc0212eSBarry Smith   c->mult_work          = 0;
13742593348eSBarry Smith   *B = C;
13753a40ed3dSBarry Smith   PetscFunctionReturn(0);
13762593348eSBarry Smith }
13772593348eSBarry Smith 
13785615d1e5SSatish Balay #undef __FUNC__
13795615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ"
138019bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
13812593348eSBarry Smith {
1382b6490206SBarry Smith   Mat_SeqBAIJ  *a;
13832593348eSBarry Smith   Mat          B;
1384de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1385b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
138635aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1387a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1388b6490206SBarry Smith   Scalar       *aa;
138919bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
13902593348eSBarry Smith 
13913a40ed3dSBarry Smith   PetscFunctionBegin;
13921a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1393de6a44a3SBarry Smith   bs2  = bs*bs;
1394b6490206SBarry Smith 
13952593348eSBarry Smith   MPI_Comm_size(comm,&size);
1396cc0fae11SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
139790ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
13980752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
1399a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
14002593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
14012593348eSBarry Smith 
1402d64ed03dSBarry Smith   if (header[3] < 0) {
1403a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1404d64ed03dSBarry Smith   }
1405d64ed03dSBarry Smith 
1406a8c6a408SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
140735aab85fSBarry Smith 
140835aab85fSBarry Smith   /*
140935aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
141035aab85fSBarry Smith     divisible by the blocksize
141135aab85fSBarry Smith   */
1412b6490206SBarry Smith   mbs        = M/bs;
141335aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
141435aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
141535aab85fSBarry Smith   else                  mbs++;
141635aab85fSBarry Smith   if (extra_rows) {
1417537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
141835aab85fSBarry Smith   }
1419b6490206SBarry Smith 
14202593348eSBarry Smith   /* read in row lengths */
142135aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
14220752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
142335aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
14242593348eSBarry Smith 
1425b6490206SBarry Smith   /* read in column indices */
142635aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
14270752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr);
142835aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1429b6490206SBarry Smith 
1430b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1431b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1432b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
143335aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
143435aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
143535aab85fSBarry Smith   masked = mask + mbs;
1436b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1437b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
143835aab85fSBarry Smith     nmask = 0;
1439b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1440b6490206SBarry Smith       kmax = rowlengths[rowcount];
1441b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
144235aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
144335aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1444b6490206SBarry Smith       }
1445b6490206SBarry Smith       rowcount++;
1446b6490206SBarry Smith     }
144735aab85fSBarry Smith     browlengths[i] += nmask;
144835aab85fSBarry Smith     /* zero out the mask elements we set */
144935aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1450b6490206SBarry Smith   }
1451b6490206SBarry Smith 
14522593348eSBarry Smith   /* create our matrix */
14533eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
14542593348eSBarry Smith   B = *A;
1455b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
14562593348eSBarry Smith 
14572593348eSBarry Smith   /* set matrix "i" values */
1458de6a44a3SBarry Smith   a->i[0] = 0;
1459b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1460b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1461b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
14622593348eSBarry Smith   }
1463b6490206SBarry Smith   a->nz         = 0;
1464b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
14652593348eSBarry Smith 
1466b6490206SBarry Smith   /* read in nonzero values */
146735aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
14680752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr);
146935aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1470b6490206SBarry Smith 
1471b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1472b6490206SBarry Smith   nzcount = 0; jcount = 0;
1473b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1474b6490206SBarry Smith     nzcountb = nzcount;
147535aab85fSBarry Smith     nmask    = 0;
1476b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1477b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1478b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
147935aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
148035aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1481b6490206SBarry Smith       }
1482b6490206SBarry Smith       rowcount++;
1483b6490206SBarry Smith     }
1484de6a44a3SBarry Smith     /* sort the masked values */
148577c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1486de6a44a3SBarry Smith 
1487b6490206SBarry Smith     /* set "j" values into matrix */
1488b6490206SBarry Smith     maskcount = 1;
148935aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
149035aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1491de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1492b6490206SBarry Smith     }
1493b6490206SBarry Smith     /* set "a" values into matrix */
1494de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1495b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1496b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1497b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1498de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
1499de6a44a3SBarry Smith         block  = mask[tmp] - 1;
1500de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
1501de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
1502b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1503b6490206SBarry Smith       }
1504b6490206SBarry Smith     }
150535aab85fSBarry Smith     /* zero out the mask elements we set */
150635aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1507b6490206SBarry Smith   }
1508a8c6a408SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1509b6490206SBarry Smith 
1510b6490206SBarry Smith   PetscFree(rowlengths);
1511b6490206SBarry Smith   PetscFree(browlengths);
1512b6490206SBarry Smith   PetscFree(aa);
1513b6490206SBarry Smith   PetscFree(jj);
1514b6490206SBarry Smith   PetscFree(mask);
1515b6490206SBarry Smith 
1516b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1517de6a44a3SBarry Smith 
15189c01be13SBarry Smith   ierr = MatView_Private(B); CHKERRQ(ierr);
15193a40ed3dSBarry Smith   PetscFunctionReturn(0);
15202593348eSBarry Smith }
15212593348eSBarry Smith 
15222593348eSBarry Smith 
15232593348eSBarry Smith 
1524