xref: /petsc/src/mat/impls/baij/seq/baij.c (revision e20fef11c1a1457fc77d0de77d14a2dd6b3370e4)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*e20fef11SSatish Balay static char vcid[] = "$Id: baij.c,v 1.137 1998/05/29 20:37:31 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"
105e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
1062d61bbb3SSatish Balay {
1072d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
108e51c0b9cSSatish Balay   int         ierr;
1092d61bbb3SSatish Balay 
11094d884c6SBarry Smith   if (--A->refct > 0) PetscFunctionReturn(0);
11194d884c6SBarry Smith 
11294d884c6SBarry Smith   if (A->mapping) {
11394d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr);
11494d884c6SBarry Smith   }
11594d884c6SBarry Smith   if (A->bmapping) {
11694d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr);
11794d884c6SBarry Smith   }
1182d61bbb3SSatish Balay #if defined(USE_PETSC_LOG)
119e1311b90SBarry Smith   PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
1202d61bbb3SSatish Balay #endif
1212d61bbb3SSatish Balay   PetscFree(a->a);
1222d61bbb3SSatish Balay   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
1232d61bbb3SSatish Balay   if (a->diag) PetscFree(a->diag);
1242d61bbb3SSatish Balay   if (a->ilen) PetscFree(a->ilen);
1252d61bbb3SSatish Balay   if (a->imax) PetscFree(a->imax);
1262d61bbb3SSatish Balay   if (a->solve_work) PetscFree(a->solve_work);
1272d61bbb3SSatish Balay   if (a->mult_work) PetscFree(a->mult_work);
128e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
1292d61bbb3SSatish Balay   PetscFree(a);
1302d61bbb3SSatish Balay   PLogObjectDestroy(A);
1312d61bbb3SSatish Balay   PetscHeaderDestroy(A);
1322d61bbb3SSatish Balay   PetscFunctionReturn(0);
1332d61bbb3SSatish Balay }
1342d61bbb3SSatish Balay 
1352d61bbb3SSatish Balay #undef __FUNC__
1362d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ"
1372d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1382d61bbb3SSatish Balay {
1392d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1402d61bbb3SSatish Balay 
1412d61bbb3SSatish Balay   PetscFunctionBegin;
1422d61bbb3SSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
1432d61bbb3SSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
1442d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
1452d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
1462d61bbb3SSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
1472d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
1482d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
1492d61bbb3SSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
1502d61bbb3SSatish Balay   else if (op == MAT_ROWS_SORTED ||
1512d61bbb3SSatish Balay            op == MAT_ROWS_UNSORTED ||
1522d61bbb3SSatish Balay            op == MAT_SYMMETRIC ||
1532d61bbb3SSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
1542d61bbb3SSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
155b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
156b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE) {
157981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
1582d61bbb3SSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
1592d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
1602d61bbb3SSatish Balay   } else {
1612d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
1622d61bbb3SSatish Balay   }
1632d61bbb3SSatish Balay   PetscFunctionReturn(0);
1642d61bbb3SSatish Balay }
1652d61bbb3SSatish Balay 
1662d61bbb3SSatish Balay 
1672d61bbb3SSatish Balay #undef __FUNC__
1682d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ"
1692d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
1702d61bbb3SSatish Balay {
1712d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1722d61bbb3SSatish Balay 
1732d61bbb3SSatish Balay   PetscFunctionBegin;
1742d61bbb3SSatish Balay   if (m) *m = a->m;
1752d61bbb3SSatish Balay   if (n) *n = a->n;
1762d61bbb3SSatish Balay   PetscFunctionReturn(0);
1772d61bbb3SSatish Balay }
1782d61bbb3SSatish Balay 
1792d61bbb3SSatish Balay #undef __FUNC__
1802d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ"
1812d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
1822d61bbb3SSatish Balay {
1832d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1842d61bbb3SSatish Balay 
1852d61bbb3SSatish Balay   PetscFunctionBegin;
1862d61bbb3SSatish Balay   *m = 0; *n = a->m;
1872d61bbb3SSatish Balay   PetscFunctionReturn(0);
1882d61bbb3SSatish Balay }
1892d61bbb3SSatish Balay 
1902d61bbb3SSatish Balay #undef __FUNC__
1912d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ"
1922d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
1932d61bbb3SSatish Balay {
1942d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1952d61bbb3SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
1962d61bbb3SSatish Balay   Scalar      *aa,*v_i,*aa_i;
1972d61bbb3SSatish Balay 
1982d61bbb3SSatish Balay   PetscFunctionBegin;
1992d61bbb3SSatish Balay   bs  = a->bs;
2002d61bbb3SSatish Balay   ai  = a->i;
2012d61bbb3SSatish Balay   aj  = a->j;
2022d61bbb3SSatish Balay   aa  = a->a;
2032d61bbb3SSatish Balay   bs2 = a->bs2;
2042d61bbb3SSatish Balay 
2052d61bbb3SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
2062d61bbb3SSatish Balay 
2072d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2082d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2092d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2102d61bbb3SSatish Balay   *nz = bs*M;
2112d61bbb3SSatish Balay 
2122d61bbb3SSatish Balay   if (v) {
2132d61bbb3SSatish Balay     *v = 0;
2142d61bbb3SSatish Balay     if (*nz) {
2152d61bbb3SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
2162d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2172d61bbb3SSatish Balay         v_i  = *v + i*bs;
2182d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2192d61bbb3SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
2202d61bbb3SSatish Balay       }
2212d61bbb3SSatish Balay     }
2222d61bbb3SSatish Balay   }
2232d61bbb3SSatish Balay 
2242d61bbb3SSatish Balay   if (idx) {
2252d61bbb3SSatish Balay     *idx = 0;
2262d61bbb3SSatish Balay     if (*nz) {
2272d61bbb3SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
2282d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2292d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2302d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2312d61bbb3SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
2322d61bbb3SSatish Balay       }
2332d61bbb3SSatish Balay     }
2342d61bbb3SSatish Balay   }
2352d61bbb3SSatish Balay   PetscFunctionReturn(0);
2362d61bbb3SSatish Balay }
2372d61bbb3SSatish Balay 
2382d61bbb3SSatish Balay #undef __FUNC__
2392d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ"
2402d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2412d61bbb3SSatish Balay {
2422d61bbb3SSatish Balay   PetscFunctionBegin;
2432d61bbb3SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
2442d61bbb3SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
2452d61bbb3SSatish Balay   PetscFunctionReturn(0);
2462d61bbb3SSatish Balay }
2472d61bbb3SSatish Balay 
2482d61bbb3SSatish Balay #undef __FUNC__
2492d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ"
2502d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
2512d61bbb3SSatish Balay {
2522d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
2532d61bbb3SSatish Balay   Mat         C;
2542d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
2552d61bbb3SSatish Balay   int         *rows,*cols,bs2=a->bs2;
2562d61bbb3SSatish Balay   Scalar      *array=a->a;
2572d61bbb3SSatish Balay 
2582d61bbb3SSatish Balay   PetscFunctionBegin;
2592d61bbb3SSatish Balay   if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
2602d61bbb3SSatish Balay   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
2612d61bbb3SSatish Balay   PetscMemzero(col,(1+nbs)*sizeof(int));
2622d61bbb3SSatish Balay 
2632d61bbb3SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
2642d61bbb3SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
2652d61bbb3SSatish Balay   PetscFree(col);
2662d61bbb3SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
2672d61bbb3SSatish Balay   cols = rows + bs;
2682d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2692d61bbb3SSatish Balay     cols[0] = i*bs;
2702d61bbb3SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
2712d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
2722d61bbb3SSatish Balay     for ( j=0; j<len; j++ ) {
2732d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
2742d61bbb3SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
2752d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
2762d61bbb3SSatish Balay       array += bs2;
2772d61bbb3SSatish Balay     }
2782d61bbb3SSatish Balay   }
2792d61bbb3SSatish Balay   PetscFree(rows);
2802d61bbb3SSatish Balay 
2812d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2822d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2832d61bbb3SSatish Balay 
2842d61bbb3SSatish Balay   if (B != PETSC_NULL) {
2852d61bbb3SSatish Balay     *B = C;
2862d61bbb3SSatish Balay   } else {
287f830108cSBarry Smith     PetscOps       *Abops;
288f830108cSBarry Smith     struct _MatOps *Aops;
289f830108cSBarry Smith 
2902d61bbb3SSatish Balay     /* This isn't really an in-place transpose */
2912d61bbb3SSatish Balay     PetscFree(a->a);
2922d61bbb3SSatish Balay     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
2932d61bbb3SSatish Balay     if (a->diag) PetscFree(a->diag);
2942d61bbb3SSatish Balay     if (a->ilen) PetscFree(a->ilen);
2952d61bbb3SSatish Balay     if (a->imax) PetscFree(a->imax);
2962d61bbb3SSatish Balay     if (a->solve_work) PetscFree(a->solve_work);
2972d61bbb3SSatish Balay     PetscFree(a);
298f830108cSBarry Smith 
299f830108cSBarry Smith     /*
300f830108cSBarry Smith        This is horrible, horrible code. We need to keep the
301f830108cSBarry Smith       A pointers for the bops and ops but copy everything
302f830108cSBarry Smith       else from C.
303f830108cSBarry Smith     */
304f830108cSBarry Smith     Abops    = A->bops;
305f830108cSBarry Smith     Aops     = A->ops;
3062d61bbb3SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
307f830108cSBarry Smith     A->bops  = Abops;
308f830108cSBarry Smith     A->ops   = Aops;
30927a8da17SBarry Smith     A->qlist = 0;
310f830108cSBarry Smith 
3112d61bbb3SSatish Balay     PetscHeaderDestroy(C);
3122d61bbb3SSatish Balay   }
3132d61bbb3SSatish Balay   PetscFunctionReturn(0);
3142d61bbb3SSatish Balay }
3152d61bbb3SSatish Balay 
3162d61bbb3SSatish Balay 
3172d61bbb3SSatish Balay 
3183b2fbd54SBarry Smith 
3195615d1e5SSatish Balay #undef __FUNC__
320d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary"
321b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
3222593348eSBarry Smith {
323b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3249df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
325b6490206SBarry Smith   Scalar      *aa;
3262593348eSBarry Smith 
3273a40ed3dSBarry Smith   PetscFunctionBegin;
32890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
3292593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
3302593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3313b2fbd54SBarry Smith 
3322593348eSBarry Smith   col_lens[1] = a->m;
3332593348eSBarry Smith   col_lens[2] = a->n;
3347e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3352593348eSBarry Smith 
3362593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
337b6490206SBarry Smith   count = 0;
338b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
339b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
340b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
341b6490206SBarry Smith     }
3422593348eSBarry Smith   }
3430752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3442593348eSBarry Smith   PetscFree(col_lens);
3452593348eSBarry Smith 
3462593348eSBarry Smith   /* store column indices (zero start index) */
34766963ce1SSatish Balay   jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj);
348b6490206SBarry Smith   count = 0;
349b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
350b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
351b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
352b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
353b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3542593348eSBarry Smith         }
3552593348eSBarry Smith       }
356b6490206SBarry Smith     }
357b6490206SBarry Smith   }
3580752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr);
359b6490206SBarry Smith   PetscFree(jj);
3602593348eSBarry Smith 
3612593348eSBarry Smith   /* store nonzero values */
36266963ce1SSatish Balay   aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa);
363b6490206SBarry Smith   count = 0;
364b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
365b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
366b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
367b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
3687e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
369b6490206SBarry Smith         }
370b6490206SBarry Smith       }
371b6490206SBarry Smith     }
372b6490206SBarry Smith   }
3730752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
374b6490206SBarry Smith   PetscFree(aa);
3753a40ed3dSBarry Smith   PetscFunctionReturn(0);
3762593348eSBarry Smith }
3772593348eSBarry Smith 
3785615d1e5SSatish Balay #undef __FUNC__
379d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII"
380b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
3812593348eSBarry Smith {
382b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3839df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
3842593348eSBarry Smith   FILE        *fd;
3852593348eSBarry Smith   char        *outputname;
3862593348eSBarry Smith 
3873a40ed3dSBarry Smith   PetscFunctionBegin;
38890ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
3892593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
39090ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
391639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
3927ddc982cSLois Curfman McInnes     fprintf(fd,"  block size is %d\n",bs);
3932593348eSBarry Smith   }
394639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
395a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
3962593348eSBarry Smith   }
397639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
39844cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
39944cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
40044cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
40144cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
40244cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
4033a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
404*e20fef11SSatish Balay           if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0)
40544cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
406*e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
407*e20fef11SSatish Balay           else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0)
408766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
409*e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));
410*e20fef11SSatish Balay           else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0)
411*e20fef11SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
41244cd7ae7SLois Curfman McInnes #else
41344cd7ae7SLois Curfman McInnes           if (a->a[bs2*k + l*bs + j] != 0.0)
41444cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
41544cd7ae7SLois Curfman McInnes #endif
41644cd7ae7SLois Curfman McInnes           }
41744cd7ae7SLois Curfman McInnes         }
41844cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
41944cd7ae7SLois Curfman McInnes       }
42044cd7ae7SLois Curfman McInnes     }
42144cd7ae7SLois Curfman McInnes   }
4222593348eSBarry Smith   else {
423b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
424b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
425b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
426b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
427b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
4283a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
429*e20fef11SSatish Balay           if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) {
43088685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
431*e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
43288685aaeSLois Curfman McInnes           }
433*e20fef11SSatish Balay           else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) {
434766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
435*e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));
436766eeae4SLois Curfman McInnes           }
43788685aaeSLois Curfman McInnes           else {
438*e20fef11SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
43988685aaeSLois Curfman McInnes           }
44088685aaeSLois Curfman McInnes #else
4417e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
44288685aaeSLois Curfman McInnes #endif
4432593348eSBarry Smith           }
4442593348eSBarry Smith         }
4452593348eSBarry Smith         fprintf(fd,"\n");
4462593348eSBarry Smith       }
4472593348eSBarry Smith     }
448b6490206SBarry Smith   }
4492593348eSBarry Smith   fflush(fd);
4503a40ed3dSBarry Smith   PetscFunctionReturn(0);
4512593348eSBarry Smith }
4522593348eSBarry Smith 
4535615d1e5SSatish Balay #undef __FUNC__
454d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw"
4553270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
4563270192aSSatish Balay {
4573270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
4583270192aSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2;
4593270192aSSatish Balay   double       xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
4603270192aSSatish Balay   Scalar       *aa;
4613270192aSSatish Balay   Draw         draw;
4623270192aSSatish Balay   DrawButton   button;
4633270192aSSatish Balay   PetscTruth   isnull;
4643270192aSSatish Balay 
4653a40ed3dSBarry Smith   PetscFunctionBegin;
4663a40ed3dSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr);
4673a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
4683270192aSSatish Balay 
4693270192aSSatish Balay   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
4703270192aSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
4713270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
4723270192aSSatish Balay   /* loop over matrix elements drawing boxes */
4733270192aSSatish Balay   color = DRAW_BLUE;
4743270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4753270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4763270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4773270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4783270192aSSatish Balay       aa = a->a + j*bs2;
4793270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4803270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4813270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
482b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4833270192aSSatish Balay         }
4843270192aSSatish Balay       }
4853270192aSSatish Balay     }
4863270192aSSatish Balay   }
4873270192aSSatish Balay   color = DRAW_CYAN;
4883270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4893270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4903270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4913270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4923270192aSSatish Balay       aa = a->a + j*bs2;
4933270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4943270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4953270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
496b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4973270192aSSatish Balay         }
4983270192aSSatish Balay       }
4993270192aSSatish Balay     }
5003270192aSSatish Balay   }
5013270192aSSatish Balay 
5023270192aSSatish Balay   color = DRAW_RED;
5033270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5043270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5053270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5063270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5073270192aSSatish Balay       aa = a->a + j*bs2;
5083270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5093270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5103270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
511b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5123270192aSSatish Balay         }
5133270192aSSatish Balay       }
5143270192aSSatish Balay     }
5153270192aSSatish Balay   }
5163270192aSSatish Balay 
51755843e3eSBarry Smith   DrawSynchronizedFlush(draw);
5183270192aSSatish Balay   DrawGetPause(draw,&pause);
5193a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
5203270192aSSatish Balay 
5213270192aSSatish Balay   /* allow the matrix to zoom or shrink */
52255843e3eSBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
5233270192aSSatish Balay   while (button != BUTTON_RIGHT) {
52455843e3eSBarry Smith     DrawSynchronizedClear(draw);
5253270192aSSatish Balay     if (button == BUTTON_LEFT) scale = .5;
5263270192aSSatish Balay     else if (button == BUTTON_CENTER) scale = 2.;
5273270192aSSatish Balay     xl = scale*(xl + w - xc) + xc - w*scale;
5283270192aSSatish Balay     xr = scale*(xr - w - xc) + xc + w*scale;
5293270192aSSatish Balay     yl = scale*(yl + h - yc) + yc - h*scale;
5303270192aSSatish Balay     yr = scale*(yr - h - yc) + yc + h*scale;
5313270192aSSatish Balay     w *= scale; h *= scale;
5323270192aSSatish Balay     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5333270192aSSatish Balay     color = DRAW_BLUE;
5343270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5353270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5363270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5373270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5383270192aSSatish Balay         aa = a->a + j*bs2;
5393270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5403270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5413270192aSSatish Balay             if (PetscReal(*aa++) >=  0.) continue;
542b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5433270192aSSatish Balay           }
5443270192aSSatish Balay         }
5453270192aSSatish Balay       }
5463270192aSSatish Balay     }
5473270192aSSatish Balay     color = DRAW_CYAN;
5483270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5493270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5503270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5513270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5523270192aSSatish Balay         aa = a->a + j*bs2;
5533270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5543270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5553270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
556b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5573270192aSSatish Balay           }
5583270192aSSatish Balay         }
5593270192aSSatish Balay       }
5603270192aSSatish Balay     }
5613270192aSSatish Balay 
5623270192aSSatish Balay     color = DRAW_RED;
5633270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5643270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5653270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5663270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5673270192aSSatish Balay         aa = a->a + j*bs2;
5683270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5693270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5703270192aSSatish Balay             if (PetscReal(*aa++) <= 0.) continue;
571b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5723270192aSSatish Balay           }
5733270192aSSatish Balay         }
5743270192aSSatish Balay       }
5753270192aSSatish Balay     }
57655843e3eSBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
5773270192aSSatish Balay   }
5783a40ed3dSBarry Smith   PetscFunctionReturn(0);
5793270192aSSatish Balay }
5803270192aSSatish Balay 
5815615d1e5SSatish Balay #undef __FUNC__
582d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ"
583e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer)
5842593348eSBarry Smith {
58519bcc07fSBarry Smith   ViewerType  vtype;
58619bcc07fSBarry Smith   int         ierr;
5872593348eSBarry Smith 
5883a40ed3dSBarry Smith   PetscFunctionBegin;
58919bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
59019bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
591a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported");
5923a40ed3dSBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
5933a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
5943a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
5953a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
5963a40ed3dSBarry Smith   } else if (vtype == DRAW_VIEWER) {
5973a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
5985cd90555SBarry Smith   } else {
5995cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
6002593348eSBarry Smith   }
6013a40ed3dSBarry Smith   PetscFunctionReturn(0);
6022593348eSBarry Smith }
603b6490206SBarry Smith 
604cd0e1443SSatish Balay 
6055615d1e5SSatish Balay #undef __FUNC__
6062d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ"
6072d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
608cd0e1443SSatish Balay {
609cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6102d61bbb3SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
6112d61bbb3SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
6122d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6132d61bbb3SSatish Balay   Scalar     *ap, *aa = a->a, zero = 0.0;
614cd0e1443SSatish Balay 
6153a40ed3dSBarry Smith   PetscFunctionBegin;
6162d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
617cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
618a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
619a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
6202d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6212c3acbe9SBarry Smith     nrow = ailen[brow];
6222d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
623a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
624a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6252d61bbb3SSatish Balay       col  = in[l] ;
6262d61bbb3SSatish Balay       bcol = col/bs;
6272d61bbb3SSatish Balay       cidx = col%bs;
6282d61bbb3SSatish Balay       ridx = row%bs;
6292d61bbb3SSatish Balay       high = nrow;
6302d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6312d61bbb3SSatish Balay       while (high-low > 5) {
632cd0e1443SSatish Balay         t = (low+high)/2;
633cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
634cd0e1443SSatish Balay         else             low  = t;
635cd0e1443SSatish Balay       }
636cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
637cd0e1443SSatish Balay         if (rp[i] > bcol) break;
638cd0e1443SSatish Balay         if (rp[i] == bcol) {
6392d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6402d61bbb3SSatish Balay           goto finished;
641cd0e1443SSatish Balay         }
642cd0e1443SSatish Balay       }
6432d61bbb3SSatish Balay       *v++ = zero;
6442d61bbb3SSatish Balay       finished:;
645cd0e1443SSatish Balay     }
646cd0e1443SSatish Balay   }
6473a40ed3dSBarry Smith   PetscFunctionReturn(0);
648cd0e1443SSatish Balay }
649cd0e1443SSatish Balay 
6502d61bbb3SSatish Balay 
6515615d1e5SSatish Balay #undef __FUNC__
65205a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ"
65392c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
65492c4ed94SBarry Smith {
65592c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6568a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
657dd9472c6SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented;
658dd9472c6SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
6590e324ae4SSatish Balay   Scalar      *ap,*value=v,*aa=a->a,*bap;
66092c4ed94SBarry Smith 
6613a40ed3dSBarry Smith   PetscFunctionBegin;
6620e324ae4SSatish Balay   if (roworiented) {
6630e324ae4SSatish Balay     stepval = (n-1)*bs;
6640e324ae4SSatish Balay   } else {
6650e324ae4SSatish Balay     stepval = (m-1)*bs;
6660e324ae4SSatish Balay   }
66792c4ed94SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
66892c4ed94SBarry Smith     row  = im[k];
6693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
670a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
671a8c6a408SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
67292c4ed94SBarry Smith #endif
67392c4ed94SBarry Smith     rp   = aj + ai[row];
67492c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
67592c4ed94SBarry Smith     rmax = imax[row];
67692c4ed94SBarry Smith     nrow = ailen[row];
67792c4ed94SBarry Smith     low  = 0;
67892c4ed94SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
6793a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
680a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
681a8c6a408SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
68292c4ed94SBarry Smith #endif
68392c4ed94SBarry Smith       col = in[l];
68492c4ed94SBarry Smith       if (roworiented) {
6850e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
6860e324ae4SSatish Balay       } else {
6870e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
68892c4ed94SBarry Smith       }
68992c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
69092c4ed94SBarry Smith       while (high-low > 7) {
69192c4ed94SBarry Smith         t = (low+high)/2;
69292c4ed94SBarry Smith         if (rp[t] > col) high = t;
69392c4ed94SBarry Smith         else             low  = t;
69492c4ed94SBarry Smith       }
69592c4ed94SBarry Smith       for ( i=low; i<high; i++ ) {
69692c4ed94SBarry Smith         if (rp[i] > col) break;
69792c4ed94SBarry Smith         if (rp[i] == col) {
6988a84c255SSatish Balay           bap  = ap +  bs2*i;
6990e324ae4SSatish Balay           if (roworiented) {
7008a84c255SSatish Balay             if (is == ADD_VALUES) {
701dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
702dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7038a84c255SSatish Balay                   bap[jj] += *value++;
704dd9472c6SBarry Smith                 }
705dd9472c6SBarry Smith               }
7060e324ae4SSatish Balay             } else {
707dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
708dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7090e324ae4SSatish Balay                   bap[jj] = *value++;
7108a84c255SSatish Balay                 }
711dd9472c6SBarry Smith               }
712dd9472c6SBarry Smith             }
7130e324ae4SSatish Balay           } else {
7140e324ae4SSatish Balay             if (is == ADD_VALUES) {
715dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
716dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7170e324ae4SSatish Balay                   *bap++ += *value++;
718dd9472c6SBarry Smith                 }
719dd9472c6SBarry Smith               }
7200e324ae4SSatish Balay             } else {
721dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
722dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7230e324ae4SSatish Balay                   *bap++  = *value++;
7240e324ae4SSatish Balay                 }
7258a84c255SSatish Balay               }
726dd9472c6SBarry Smith             }
727dd9472c6SBarry Smith           }
728f1241b54SBarry Smith           goto noinsert2;
72992c4ed94SBarry Smith         }
73092c4ed94SBarry Smith       }
73189280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
732a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
73392c4ed94SBarry Smith       if (nrow >= rmax) {
73492c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
73592c4ed94SBarry Smith         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
73692c4ed94SBarry Smith         Scalar *new_a;
73792c4ed94SBarry Smith 
738a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
73989280ab3SLois Curfman McInnes 
74092c4ed94SBarry Smith         /* malloc new storage space */
74192c4ed94SBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
74292c4ed94SBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
74392c4ed94SBarry Smith         new_j   = (int *) (new_a + bs2*new_nz);
74492c4ed94SBarry Smith         new_i   = new_j + new_nz;
74592c4ed94SBarry Smith 
74692c4ed94SBarry Smith         /* copy over old data into new slots */
74792c4ed94SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
74892c4ed94SBarry Smith         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
74992c4ed94SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));
75092c4ed94SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow);
751dd9472c6SBarry Smith         PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));
75292c4ed94SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(Scalar));
75392c4ed94SBarry Smith         PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
754dd9472c6SBarry Smith         PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(Scalar));
75592c4ed94SBarry Smith         /* free up old matrix storage */
75692c4ed94SBarry Smith         PetscFree(a->a);
75792c4ed94SBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
75892c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
75992c4ed94SBarry Smith         a->singlemalloc = 1;
76092c4ed94SBarry Smith 
76192c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
76292c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
76392c4ed94SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
76492c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
76592c4ed94SBarry Smith         a->reallocs++;
76692c4ed94SBarry Smith         a->nz++;
76792c4ed94SBarry Smith       }
76892c4ed94SBarry Smith       N = nrow++ - 1;
76992c4ed94SBarry Smith       /* shift up all the later entries in this row */
77092c4ed94SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
77192c4ed94SBarry Smith         rp[ii+1] = rp[ii];
77292c4ed94SBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
77392c4ed94SBarry Smith       }
77492c4ed94SBarry Smith       if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
77592c4ed94SBarry Smith       rp[i] = col;
7768a84c255SSatish Balay       bap   = ap +  bs2*i;
7770e324ae4SSatish Balay       if (roworiented) {
778dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval) {
779dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs ) {
7800e324ae4SSatish Balay             bap[jj] = *value++;
781dd9472c6SBarry Smith           }
782dd9472c6SBarry Smith         }
7830e324ae4SSatish Balay       } else {
784dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval ) {
785dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++ ) {
7860e324ae4SSatish Balay             *bap++  = *value++;
7870e324ae4SSatish Balay           }
788dd9472c6SBarry Smith         }
789dd9472c6SBarry Smith       }
790f1241b54SBarry Smith       noinsert2:;
79192c4ed94SBarry Smith       low = i;
79292c4ed94SBarry Smith     }
79392c4ed94SBarry Smith     ailen[row] = nrow;
79492c4ed94SBarry Smith   }
7953a40ed3dSBarry Smith   PetscFunctionReturn(0);
79692c4ed94SBarry Smith }
79792c4ed94SBarry Smith 
798f2501298SSatish Balay 
7995615d1e5SSatish Balay #undef __FUNC__
8005615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ"
8018f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
802584200bdSSatish Balay {
803584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
804584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
805a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
80643ee02c3SBarry Smith   int        mbs = a->mbs, bs2 = a->bs2,rmax = 0;
807584200bdSSatish Balay   Scalar     *aa = a->a, *ap;
808584200bdSSatish Balay 
8093a40ed3dSBarry Smith   PetscFunctionBegin;
8103a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
811584200bdSSatish Balay 
81243ee02c3SBarry Smith   if (m) rmax = ailen[0];
813584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
814584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
815584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
816d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
817584200bdSSatish Balay     if (fshift) {
818a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
819584200bdSSatish Balay       N = ailen[i];
820584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
821584200bdSSatish Balay         ip[j-fshift] = ip[j];
8227e67e3f9SSatish Balay         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar));
823584200bdSSatish Balay       }
824584200bdSSatish Balay     }
825584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
826584200bdSSatish Balay   }
827584200bdSSatish Balay   if (mbs) {
828584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
829584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
830584200bdSSatish Balay   }
831584200bdSSatish Balay   /* reset ilen and imax for each row */
832584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
833584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
834584200bdSSatish Balay   }
835a7c10996SSatish Balay   a->nz = ai[mbs];
836584200bdSSatish Balay 
837584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
838584200bdSSatish Balay   if (fshift && a->diag) {
839584200bdSSatish Balay     PetscFree(a->diag);
840584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
841584200bdSSatish Balay     a->diag = 0;
842584200bdSSatish Balay   }
8433d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
844219d9a1aSLois Curfman McInnes            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
8453d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
846584200bdSSatish Balay            a->reallocs);
847d402145bSBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
848e2f3b5e9SSatish Balay   a->reallocs          = 0;
8494e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
8504e220ebcSLois Curfman McInnes 
8513a40ed3dSBarry Smith   PetscFunctionReturn(0);
852584200bdSSatish Balay }
853584200bdSSatish Balay 
8545a838052SSatish Balay 
855d9b7c43dSSatish Balay /* idx should be of length atlease bs */
8565615d1e5SSatish Balay #undef __FUNC__
8575615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block"
858d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg)
859d9b7c43dSSatish Balay {
860d9b7c43dSSatish Balay   int i,row;
8613a40ed3dSBarry Smith 
8623a40ed3dSBarry Smith   PetscFunctionBegin;
863d9b7c43dSSatish Balay   row = idx[0];
8643a40ed3dSBarry Smith   if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
865d9b7c43dSSatish Balay 
866d9b7c43dSSatish Balay   for ( i=1; i<bs; i++ ) {
8673a40ed3dSBarry Smith     if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
868d9b7c43dSSatish Balay   }
869d9b7c43dSSatish Balay   *flg = PETSC_TRUE;
8703a40ed3dSBarry Smith   PetscFunctionReturn(0);
871d9b7c43dSSatish Balay }
872d9b7c43dSSatish Balay 
8735615d1e5SSatish Balay #undef __FUNC__
8745615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ"
8758f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
876d9b7c43dSSatish Balay {
877d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
878d9b7c43dSSatish Balay   IS          is_local;
879d9b7c43dSSatish Balay   int         ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2;
880d9b7c43dSSatish Balay   PetscTruth  flg;
881d9b7c43dSSatish Balay   Scalar      zero = 0.0,*aa;
882d9b7c43dSSatish Balay 
8833a40ed3dSBarry Smith   PetscFunctionBegin;
884d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
885d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
886d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
887537820f0SBarry Smith   ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr);
888d9b7c43dSSatish Balay   ierr = ISSort(is_local); CHKERRQ(ierr);
889d9b7c43dSSatish Balay   ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr);
890d9b7c43dSSatish Balay 
891d9b7c43dSSatish Balay   i = 0;
892d9b7c43dSSatish Balay   while (i < is_n) {
893a8c6a408SBarry Smith     if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
894d9b7c43dSSatish Balay     flg = PETSC_FALSE;
895d9b7c43dSSatish Balay     if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); }
896d9b7c43dSSatish Balay     count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs;
897d9b7c43dSSatish Balay     aa    = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs);
8982d61bbb3SSatish Balay     if (flg) { /* There exists a block of rows to be Zerowed */
899a07cd24cSSatish Balay       if (baij->ilen[rows[i]/bs] > 0) {
9002d61bbb3SSatish Balay         PetscMemzero(aa,count*bs*sizeof(Scalar));
9012d61bbb3SSatish Balay         baij->ilen[rows[i]/bs] = 1;
9022d61bbb3SSatish Balay         baij->j[baij->i[rows[i]/bs]] = rows[i]/bs;
903a07cd24cSSatish Balay       }
9042d61bbb3SSatish Balay       i += bs;
9052d61bbb3SSatish Balay     } else { /* Zero out only the requested row */
906d9b7c43dSSatish Balay       for ( j=0; j<count; j++ ) {
907d9b7c43dSSatish Balay         aa[0] = zero;
908d9b7c43dSSatish Balay         aa+=bs;
909d9b7c43dSSatish Balay       }
910d9b7c43dSSatish Balay       i++;
911d9b7c43dSSatish Balay     }
912d9b7c43dSSatish Balay   }
913d9b7c43dSSatish Balay   if (diag) {
914d9b7c43dSSatish Balay     for ( j=0; j<is_n; j++ ) {
915f830108cSBarry Smith       ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
9162d61bbb3SSatish Balay       /* ierr = MatSetValues(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); */
917d9b7c43dSSatish Balay     }
918d9b7c43dSSatish Balay   }
919d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr);
920d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr);
921d9b7c43dSSatish Balay   ierr = ISDestroy(is_local); CHKERRQ(ierr);
9229a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9233a40ed3dSBarry Smith   PetscFunctionReturn(0);
924d9b7c43dSSatish Balay }
9251c351548SSatish Balay 
9265615d1e5SSatish Balay #undef __FUNC__
9272d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ"
9282d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
9292d61bbb3SSatish Balay {
9302d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9312d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
9322d61bbb3SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
9332d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
9342d61bbb3SSatish Balay   int          ridx,cidx,bs2=a->bs2;
9352d61bbb3SSatish Balay   Scalar      *ap,value,*aa=a->a,*bap;
9362d61bbb3SSatish Balay 
9372d61bbb3SSatish Balay   PetscFunctionBegin;
9382d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
9392d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
9402d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9412d61bbb3SSatish Balay     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
9422d61bbb3SSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
9432d61bbb3SSatish Balay #endif
9442d61bbb3SSatish Balay     rp   = aj + ai[brow];
9452d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
9462d61bbb3SSatish Balay     rmax = imax[brow];
9472d61bbb3SSatish Balay     nrow = ailen[brow];
9482d61bbb3SSatish Balay     low  = 0;
9492d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
9502d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9512d61bbb3SSatish Balay       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
9522d61bbb3SSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
9532d61bbb3SSatish Balay #endif
9542d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
9552d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
9562d61bbb3SSatish Balay       if (roworiented) {
9572d61bbb3SSatish Balay         value = *v++;
9582d61bbb3SSatish Balay       } else {
9592d61bbb3SSatish Balay         value = v[k + l*m];
9602d61bbb3SSatish Balay       }
9612d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
9622d61bbb3SSatish Balay       while (high-low > 7) {
9632d61bbb3SSatish Balay         t = (low+high)/2;
9642d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
9652d61bbb3SSatish Balay         else              low  = t;
9662d61bbb3SSatish Balay       }
9672d61bbb3SSatish Balay       for ( i=low; i<high; i++ ) {
9682d61bbb3SSatish Balay         if (rp[i] > bcol) break;
9692d61bbb3SSatish Balay         if (rp[i] == bcol) {
9702d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
9712d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
9722d61bbb3SSatish Balay           else                  *bap  = value;
9732d61bbb3SSatish Balay           goto noinsert1;
9742d61bbb3SSatish Balay         }
9752d61bbb3SSatish Balay       }
9762d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
9772d61bbb3SSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9782d61bbb3SSatish Balay       if (nrow >= rmax) {
9792d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
9802d61bbb3SSatish Balay         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
9812d61bbb3SSatish Balay         Scalar *new_a;
9822d61bbb3SSatish Balay 
9832d61bbb3SSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9842d61bbb3SSatish Balay 
9852d61bbb3SSatish Balay         /* Malloc new storage space */
9862d61bbb3SSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
9872d61bbb3SSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
9882d61bbb3SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
9892d61bbb3SSatish Balay         new_i   = new_j + new_nz;
9902d61bbb3SSatish Balay 
9912d61bbb3SSatish Balay         /* copy over old data into new slots */
9922d61bbb3SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
9932d61bbb3SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
9942d61bbb3SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
9952d61bbb3SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
9962d61bbb3SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
9972d61bbb3SSatish Balay                                                            len*sizeof(int));
9982d61bbb3SSatish Balay         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar));
9992d61bbb3SSatish Balay         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
10002d61bbb3SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
10012d61bbb3SSatish Balay                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar));
10022d61bbb3SSatish Balay         /* free up old matrix storage */
10032d61bbb3SSatish Balay         PetscFree(a->a);
10042d61bbb3SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
10052d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
10062d61bbb3SSatish Balay         a->singlemalloc = 1;
10072d61bbb3SSatish Balay 
10082d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10092d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
10102d61bbb3SSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
10112d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10122d61bbb3SSatish Balay         a->reallocs++;
10132d61bbb3SSatish Balay         a->nz++;
10142d61bbb3SSatish Balay       }
10152d61bbb3SSatish Balay       N = nrow++ - 1;
10162d61bbb3SSatish Balay       /* shift up all the later entries in this row */
10172d61bbb3SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
10182d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
10192d61bbb3SSatish Balay         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
10202d61bbb3SSatish Balay       }
10212d61bbb3SSatish Balay       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
10222d61bbb3SSatish Balay       rp[i]                      = bcol;
10232d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
10242d61bbb3SSatish Balay       noinsert1:;
10252d61bbb3SSatish Balay       low = i;
10262d61bbb3SSatish Balay     }
10272d61bbb3SSatish Balay     ailen[brow] = nrow;
10282d61bbb3SSatish Balay   }
10292d61bbb3SSatish Balay   PetscFunctionReturn(0);
10302d61bbb3SSatish Balay }
10312d61bbb3SSatish Balay 
10322d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
10332d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
10342d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
1035d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*);
10362d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**);
10372d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec);
10382d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
10392d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat);
10402d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *);
10412d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*);
10422d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec);
10432d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
10442d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
10452d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat);
10462d61bbb3SSatish Balay 
10472d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
10482d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
10492d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
10502d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
10512d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
10522d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
10532d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
10542d61bbb3SSatish Balay 
10552d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
10562d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
10572d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
10582d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
10592d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
10602d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
10612d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
10622d61bbb3SSatish Balay 
10632d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
10642d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
10652d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
10662d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
10672d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
10682d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
10692d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
10702d61bbb3SSatish Balay 
10712d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
10722d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
10732d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
10742d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
10752d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
10762d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
10772d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
10782d61bbb3SSatish Balay 
10792d61bbb3SSatish Balay /*
10802d61bbb3SSatish Balay      note: This can only work for identity for row and col. It would
10812d61bbb3SSatish Balay    be good to check this and otherwise generate an error.
10822d61bbb3SSatish Balay */
10832d61bbb3SSatish Balay #undef __FUNC__
10842d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ"
10852d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
10862d61bbb3SSatish Balay {
10872d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
10882d61bbb3SSatish Balay   Mat         outA;
10892d61bbb3SSatish Balay   int         ierr;
10902d61bbb3SSatish Balay 
10912d61bbb3SSatish Balay   PetscFunctionBegin;
10922d61bbb3SSatish Balay   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
10932d61bbb3SSatish Balay 
10942d61bbb3SSatish Balay   outA          = inA;
10952d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
10962d61bbb3SSatish Balay   a->row        = row;
10972d61bbb3SSatish Balay   a->col        = col;
10982d61bbb3SSatish Balay 
1099e51c0b9cSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1100e51c0b9cSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
11011e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1102e51c0b9cSSatish Balay 
11032d61bbb3SSatish Balay   if (!a->solve_work) {
11042d61bbb3SSatish Balay     a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
11052d61bbb3SSatish Balay     PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
11062d61bbb3SSatish Balay   }
11072d61bbb3SSatish Balay 
11082d61bbb3SSatish Balay   if (!a->diag) {
11092d61bbb3SSatish Balay     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
11102d61bbb3SSatish Balay   }
11112d61bbb3SSatish Balay   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
11122d61bbb3SSatish Balay 
11132d61bbb3SSatish Balay   /*
11142d61bbb3SSatish Balay       Blocksize 4 has a special faster solver for ILU(0) factorization
11152d61bbb3SSatish Balay     with natural ordering
11162d61bbb3SSatish Balay   */
11172d61bbb3SSatish Balay   if (a->bs == 4) {
1118f830108cSBarry Smith     inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering;
11192d61bbb3SSatish Balay   }
11202d61bbb3SSatish Balay 
11212d61bbb3SSatish Balay   PetscFunctionReturn(0);
11222d61bbb3SSatish Balay }
11232d61bbb3SSatish Balay #undef __FUNC__
1124d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ"
1125ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1126ba4ca20aSSatish Balay {
1127ba4ca20aSSatish Balay   static int called = 0;
1128ba4ca20aSSatish Balay   MPI_Comm   comm = A->comm;
1129ba4ca20aSSatish Balay 
11303a40ed3dSBarry Smith   PetscFunctionBegin;
11313a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
113276be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");
113376be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");
11343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1135ba4ca20aSSatish Balay }
1136d9b7c43dSSatish Balay 
113727a8da17SBarry Smith #undef __FUNC__
113827a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
113927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
114027a8da17SBarry Smith {
114127a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
114227a8da17SBarry Smith   int         i,nz,n;
114327a8da17SBarry Smith 
114427a8da17SBarry Smith   PetscFunctionBegin;
114527a8da17SBarry Smith   nz = baij->maxnz;
114627a8da17SBarry Smith   n  = baij->n;
114727a8da17SBarry Smith   for (i=0; i<nz; i++) {
114827a8da17SBarry Smith     baij->j[i] = indices[i];
114927a8da17SBarry Smith   }
115027a8da17SBarry Smith   baij->nz = nz;
115127a8da17SBarry Smith   for ( i=0; i<n; i++ ) {
115227a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
115327a8da17SBarry Smith   }
115427a8da17SBarry Smith 
115527a8da17SBarry Smith   PetscFunctionReturn(0);
115627a8da17SBarry Smith }
115727a8da17SBarry Smith 
115827a8da17SBarry Smith #undef __FUNC__
115927a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices"
116027a8da17SBarry Smith /*@
116127a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
116227a8da17SBarry Smith        in the matrix.
116327a8da17SBarry Smith 
116427a8da17SBarry Smith   Input Parameters:
116527a8da17SBarry Smith +  mat - the SeqBAIJ matrix
116627a8da17SBarry Smith -  indices - the column indices
116727a8da17SBarry Smith 
116827a8da17SBarry Smith   Notes:
116927a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
117027a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
117127a8da17SBarry Smith   of the MatSetValues() operation.
117227a8da17SBarry Smith 
117327a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
117427a8da17SBarry Smith   MatCreateSeqBAIJ().
117527a8da17SBarry Smith 
117627a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
117727a8da17SBarry Smith 
117827a8da17SBarry Smith @*/
117927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
118027a8da17SBarry Smith {
118127a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
118227a8da17SBarry Smith 
118327a8da17SBarry Smith   PetscFunctionBegin;
118427a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
118527a8da17SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);
118627a8da17SBarry Smith          CHKERRQ(ierr);
118727a8da17SBarry Smith   if (f) {
118827a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
118927a8da17SBarry Smith   } else {
119027a8da17SBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
119127a8da17SBarry Smith   }
119227a8da17SBarry Smith   PetscFunctionReturn(0);
119327a8da17SBarry Smith }
119427a8da17SBarry Smith 
11952593348eSBarry Smith /* -------------------------------------------------------------------*/
1196cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ,
11979867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
1198f44d3a6dSSatish Balay        MatMult_SeqBAIJ_N,MatMultAdd_SeqBAIJ_N,
1199faf6580fSSatish Balay        MatMultTrans_SeqBAIJ,MatMultTransAdd_SeqBAIJ,
12001a6a6d98SBarry Smith        MatSolve_SeqBAIJ_N,0,
1201b6490206SBarry Smith        0,0,
1202de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
1203b6490206SBarry Smith        0,
1204f2501298SSatish Balay        MatTranspose_SeqBAIJ,
120517e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
12069867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
1207584200bdSSatish Balay        0,MatAssemblyEnd_SeqBAIJ,
1208b6490206SBarry Smith        0,
1209d9b7c43dSSatish Balay        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,MatZeroRows_SeqBAIJ,
12107fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
1211b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
1212de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
1213d402145bSBarry Smith        0,0,
1214b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
1215b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
1216af015674SSatish Balay        MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ,
12177e67e3f9SSatish Balay        MatGetValues_SeqBAIJ,0,
1218ba4ca20aSSatish Balay        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
12193b2fbd54SBarry Smith        0,0,0,MatGetBlockSize_SeqBAIJ,
12203b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
122192c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
122292c4ed94SBarry Smith        0,0,0,0,0,0,
1223d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1224d3825aa8SBarry Smith        MatGetSubMatrix_SeqBAIJ};
12252593348eSBarry Smith 
12265615d1e5SSatish Balay #undef __FUNC__
12275615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ"
12282593348eSBarry Smith /*@C
122944cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
123044cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
123144cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
12322bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
12332bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
12342593348eSBarry Smith 
1235db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1236db81eaa0SLois Curfman McInnes 
12372593348eSBarry Smith    Input Parameters:
1238db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
1239b6490206SBarry Smith .  bs - size of block
12402593348eSBarry Smith .  m - number of rows
12412593348eSBarry Smith .  n - number of columns
1242b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1243db81eaa0SLois Curfman McInnes -  nzz - array containing the number of block nonzeros in the various block rows
12442bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
12452593348eSBarry Smith 
12462593348eSBarry Smith    Output Parameter:
12472593348eSBarry Smith .  A - the matrix
12482593348eSBarry Smith 
12490513a670SBarry Smith    Options Database Keys:
1250db81eaa0SLois Curfman McInnes .   -mat_no_unroll - uses code that does not unroll the loops in the
1251db81eaa0SLois Curfman McInnes                      block calculations (much slower)
1252db81eaa0SLois Curfman McInnes .    -mat_block_size - size of the blocks to use
12530513a670SBarry Smith 
12542593348eSBarry Smith    Notes:
125544cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
12562593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
125744cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
12582593348eSBarry Smith 
12592593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
12602593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
12612593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
12626da5968aSLois Curfman McInnes    matrices.
12632593348eSBarry Smith 
1264db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
12652593348eSBarry Smith @*/
1266b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
12672593348eSBarry Smith {
12682593348eSBarry Smith   Mat         B;
1269b6490206SBarry Smith   Mat_SeqBAIJ *b;
12703b2fbd54SBarry Smith   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size;
12713b2fbd54SBarry Smith 
12723a40ed3dSBarry Smith   PetscFunctionBegin;
12733b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1274a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1275b6490206SBarry Smith 
12760513a670SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
12770513a670SBarry Smith 
12783a40ed3dSBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1279a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
12803a40ed3dSBarry Smith   }
12812593348eSBarry Smith 
12822593348eSBarry Smith   *A = 0;
1283f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView);
12842593348eSBarry Smith   PLogObjectCreate(B);
1285b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
128644cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
1287f830108cSBarry Smith   PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps));
12881a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
12891a6a6d98SBarry Smith   if (!flg) {
12907fc0212eSBarry Smith     switch (bs) {
12917fc0212eSBarry Smith     case 1:
1292f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1293f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1294f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1295f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
12967fc0212eSBarry Smith       break;
12974eeb42bcSBarry Smith     case 2:
1298f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1299f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1300f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1301f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
13026d84be18SBarry Smith       break;
1303f327dd97SBarry Smith     case 3:
1304f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1305f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1306f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1307f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
13084eeb42bcSBarry Smith       break;
13096d84be18SBarry Smith     case 4:
1310f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1311f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1312f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1313f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
13146d84be18SBarry Smith       break;
13156d84be18SBarry Smith     case 5:
1316f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1317f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1318f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1319f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
13206d84be18SBarry Smith       break;
132148e9ddb2SSatish Balay     case 7:
1322f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1323f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1324f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
132548e9ddb2SSatish Balay       break;
13267fc0212eSBarry Smith     }
13271a6a6d98SBarry Smith   }
1328e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_SeqBAIJ;
1329e1311b90SBarry Smith   B->ops->view             = MatView_SeqBAIJ;
13302593348eSBarry Smith   B->factor           = 0;
13312593348eSBarry Smith   B->lupivotthreshold = 1.0;
133290f02eecSBarry Smith   B->mapping          = 0;
13332593348eSBarry Smith   b->row              = 0;
13342593348eSBarry Smith   b->col              = 0;
1335e51c0b9cSSatish Balay   b->icol             = 0;
13362593348eSBarry Smith   b->reallocs         = 0;
13372593348eSBarry Smith 
133844cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
133944cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1340b6490206SBarry Smith   b->mbs     = mbs;
1341f2501298SSatish Balay   b->nbs     = nbs;
1342b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
13432593348eSBarry Smith   if (nnz == PETSC_NULL) {
1344119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
13452593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1346b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1347b6490206SBarry Smith     nz = nz*mbs;
13483a40ed3dSBarry Smith   } else {
13492593348eSBarry Smith     nz = 0;
1350b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
13512593348eSBarry Smith   }
13522593348eSBarry Smith 
13532593348eSBarry Smith   /* allocate the matrix space */
13547e67e3f9SSatish Balay   len   = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int);
13552593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
13567e67e3f9SSatish Balay   PetscMemzero(b->a,nz*bs2*sizeof(Scalar));
13577e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
13582593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
13592593348eSBarry Smith   b->i  = b->j + nz;
13602593348eSBarry Smith   b->singlemalloc = 1;
13612593348eSBarry Smith 
1362de6a44a3SBarry Smith   b->i[0] = 0;
1363b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
13642593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
13652593348eSBarry Smith   }
13662593348eSBarry Smith 
1367b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1368b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1369f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1370b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
13712593348eSBarry Smith 
1372b6490206SBarry Smith   b->bs               = bs;
13739df24120SSatish Balay   b->bs2              = bs2;
1374b6490206SBarry Smith   b->mbs              = mbs;
13752593348eSBarry Smith   b->nz               = 0;
137618e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
13772593348eSBarry Smith   b->sorted           = 0;
13782593348eSBarry Smith   b->roworiented      = 1;
13792593348eSBarry Smith   b->nonew            = 0;
13802593348eSBarry Smith   b->diag             = 0;
13812593348eSBarry Smith   b->solve_work       = 0;
1382de6a44a3SBarry Smith   b->mult_work        = 0;
13832593348eSBarry Smith   b->spptr            = 0;
13844e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
13854e220ebcSLois Curfman McInnes 
13862593348eSBarry Smith   *A = B;
13872593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
13882593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
138927a8da17SBarry Smith 
139027a8da17SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
139127a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
139227a8da17SBarry Smith                                      (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
139327a8da17SBarry Smith 
13943a40ed3dSBarry Smith   PetscFunctionReturn(0);
13952593348eSBarry Smith }
13962593348eSBarry Smith 
13975615d1e5SSatish Balay #undef __FUNC__
13985615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ"
1399b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
14002593348eSBarry Smith {
14012593348eSBarry Smith   Mat         C;
1402b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
140327a8da17SBarry Smith   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1404de6a44a3SBarry Smith 
14053a40ed3dSBarry Smith   PetscFunctionBegin;
1406a8c6a408SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
14072593348eSBarry Smith 
14082593348eSBarry Smith   *B = 0;
1409f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView);
14102593348eSBarry Smith   PLogObjectCreate(C);
1411b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
1412c3c66ee5SSatish Balay   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1413e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqBAIJ;
1414e1311b90SBarry Smith   C->ops->view       = MatView_SeqBAIJ;
14152593348eSBarry Smith   C->factor     = A->factor;
14162593348eSBarry Smith   c->row        = 0;
14172593348eSBarry Smith   c->col        = 0;
14182593348eSBarry Smith   C->assembled  = PETSC_TRUE;
14192593348eSBarry Smith 
142044cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
142144cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
142244cd7ae7SLois Curfman McInnes   C->M          = a->m;
142344cd7ae7SLois Curfman McInnes   C->N          = a->n;
142444cd7ae7SLois Curfman McInnes 
1425b6490206SBarry Smith   c->bs         = a->bs;
14269df24120SSatish Balay   c->bs2        = a->bs2;
1427b6490206SBarry Smith   c->mbs        = a->mbs;
1428b6490206SBarry Smith   c->nbs        = a->nbs;
14292593348eSBarry Smith 
1430b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1431b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1432b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
14332593348eSBarry Smith     c->imax[i] = a->imax[i];
14342593348eSBarry Smith     c->ilen[i] = a->ilen[i];
14352593348eSBarry Smith   }
14362593348eSBarry Smith 
14372593348eSBarry Smith   /* allocate the matrix space */
14382593348eSBarry Smith   c->singlemalloc = 1;
14397e67e3f9SSatish Balay   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int));
14402593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
14417e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1442de6a44a3SBarry Smith   c->i  = c->j + nz;
1443b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1444b6490206SBarry Smith   if (mbs > 0) {
1445de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
14462593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
14477e67e3f9SSatish Balay       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar));
14482593348eSBarry Smith     }
14492593348eSBarry Smith   }
14502593348eSBarry Smith 
1451f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
14522593348eSBarry Smith   c->sorted      = a->sorted;
14532593348eSBarry Smith   c->roworiented = a->roworiented;
14542593348eSBarry Smith   c->nonew       = a->nonew;
14552593348eSBarry Smith 
14562593348eSBarry Smith   if (a->diag) {
1457b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1458b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1459b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
14602593348eSBarry Smith       c->diag[i] = a->diag[i];
14612593348eSBarry Smith     }
14622593348eSBarry Smith   }
14632593348eSBarry Smith   else c->diag          = 0;
14642593348eSBarry Smith   c->nz                 = a->nz;
14652593348eSBarry Smith   c->maxnz              = a->maxnz;
14662593348eSBarry Smith   c->solve_work         = 0;
14672593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
14687fc0212eSBarry Smith   c->mult_work          = 0;
14692593348eSBarry Smith   *B = C;
147027a8da17SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqBAIJSetColumnIndices_C",
147127a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
147227a8da17SBarry Smith                                      (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
14733a40ed3dSBarry Smith   PetscFunctionReturn(0);
14742593348eSBarry Smith }
14752593348eSBarry Smith 
14765615d1e5SSatish Balay #undef __FUNC__
14775615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ"
147819bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
14792593348eSBarry Smith {
1480b6490206SBarry Smith   Mat_SeqBAIJ  *a;
14812593348eSBarry Smith   Mat          B;
1482de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1483b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
148435aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1485a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1486b6490206SBarry Smith   Scalar       *aa;
148719bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
14882593348eSBarry Smith 
14893a40ed3dSBarry Smith   PetscFunctionBegin;
14901a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1491de6a44a3SBarry Smith   bs2  = bs*bs;
1492b6490206SBarry Smith 
14932593348eSBarry Smith   MPI_Comm_size(comm,&size);
1494cc0fae11SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
149590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
14960752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
1497a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
14982593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
14992593348eSBarry Smith 
1500d64ed03dSBarry Smith   if (header[3] < 0) {
1501a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1502d64ed03dSBarry Smith   }
1503d64ed03dSBarry Smith 
1504a8c6a408SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
150535aab85fSBarry Smith 
150635aab85fSBarry Smith   /*
150735aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
150835aab85fSBarry Smith     divisible by the blocksize
150935aab85fSBarry Smith   */
1510b6490206SBarry Smith   mbs        = M/bs;
151135aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
151235aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
151335aab85fSBarry Smith   else                  mbs++;
151435aab85fSBarry Smith   if (extra_rows) {
1515537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
151635aab85fSBarry Smith   }
1517b6490206SBarry Smith 
15182593348eSBarry Smith   /* read in row lengths */
151935aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
15200752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
152135aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
15222593348eSBarry Smith 
1523b6490206SBarry Smith   /* read in column indices */
152435aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
15250752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr);
152635aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1527b6490206SBarry Smith 
1528b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1529b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1530b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
153135aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
153235aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
153335aab85fSBarry Smith   masked = mask + mbs;
1534b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1535b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
153635aab85fSBarry Smith     nmask = 0;
1537b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1538b6490206SBarry Smith       kmax = rowlengths[rowcount];
1539b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
154035aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
154135aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1542b6490206SBarry Smith       }
1543b6490206SBarry Smith       rowcount++;
1544b6490206SBarry Smith     }
154535aab85fSBarry Smith     browlengths[i] += nmask;
154635aab85fSBarry Smith     /* zero out the mask elements we set */
154735aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1548b6490206SBarry Smith   }
1549b6490206SBarry Smith 
15502593348eSBarry Smith   /* create our matrix */
15513eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
15522593348eSBarry Smith   B = *A;
1553b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
15542593348eSBarry Smith 
15552593348eSBarry Smith   /* set matrix "i" values */
1556de6a44a3SBarry Smith   a->i[0] = 0;
1557b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1558b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1559b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
15602593348eSBarry Smith   }
1561b6490206SBarry Smith   a->nz         = 0;
1562b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
15632593348eSBarry Smith 
1564b6490206SBarry Smith   /* read in nonzero values */
156535aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
15660752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr);
156735aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1568b6490206SBarry Smith 
1569b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1570b6490206SBarry Smith   nzcount = 0; jcount = 0;
1571b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1572b6490206SBarry Smith     nzcountb = nzcount;
157335aab85fSBarry Smith     nmask    = 0;
1574b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1575b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1576b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
157735aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
157835aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1579b6490206SBarry Smith       }
1580b6490206SBarry Smith       rowcount++;
1581b6490206SBarry Smith     }
1582de6a44a3SBarry Smith     /* sort the masked values */
158377c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1584de6a44a3SBarry Smith 
1585b6490206SBarry Smith     /* set "j" values into matrix */
1586b6490206SBarry Smith     maskcount = 1;
158735aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
158835aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1589de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1590b6490206SBarry Smith     }
1591b6490206SBarry Smith     /* set "a" values into matrix */
1592de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1593b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1594b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1595b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1596de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
1597de6a44a3SBarry Smith         block  = mask[tmp] - 1;
1598de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
1599de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
1600b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1601b6490206SBarry Smith       }
1602b6490206SBarry Smith     }
160335aab85fSBarry Smith     /* zero out the mask elements we set */
160435aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1605b6490206SBarry Smith   }
1606a8c6a408SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1607b6490206SBarry Smith 
1608b6490206SBarry Smith   PetscFree(rowlengths);
1609b6490206SBarry Smith   PetscFree(browlengths);
1610b6490206SBarry Smith   PetscFree(aa);
1611b6490206SBarry Smith   PetscFree(jj);
1612b6490206SBarry Smith   PetscFree(mask);
1613b6490206SBarry Smith 
1614b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1615de6a44a3SBarry Smith 
16169c01be13SBarry Smith   ierr = MatView_Private(B); CHKERRQ(ierr);
16173a40ed3dSBarry Smith   PetscFunctionReturn(0);
16182593348eSBarry Smith }
16192593348eSBarry Smith 
16202593348eSBarry Smith 
16212593348eSBarry Smith 
1622