xref: /petsc/src/mat/impls/baij/seq/baij.c (revision ecc615b2352880ebb2ffd7a01b13a1b2064c02f5)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*ecc615b2SBarry Smith static char vcid[] = "$Id: baij.c,v 1.143 1998/07/16 16:00:17 bsmith Exp bsmith $";
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++ ) {
30*ecc615b2SBarry Smith     diag[i] = a->i[i+1];
31de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
32de6a44a3SBarry Smith       if (a->j[j] == i) {
33de6a44a3SBarry Smith         diag[i] = j;
34de6a44a3SBarry Smith         break;
35de6a44a3SBarry Smith       }
36de6a44a3SBarry Smith     }
37de6a44a3SBarry Smith   }
38de6a44a3SBarry Smith   a->diag = diag;
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
40de6a44a3SBarry Smith }
412593348eSBarry Smith 
422593348eSBarry Smith 
433b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
443b2fbd54SBarry Smith 
455615d1e5SSatish Balay #undef __FUNC__
465615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ"
473b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
483b2fbd54SBarry Smith                             PetscTruth *done)
493b2fbd54SBarry Smith {
503b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
513b2fbd54SBarry Smith   int         ierr, n = a->mbs,i;
523b2fbd54SBarry Smith 
533a40ed3dSBarry Smith   PetscFunctionBegin;
543b2fbd54SBarry Smith   *nn = n;
553a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
563b2fbd54SBarry Smith   if (symmetric) {
573b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
583b2fbd54SBarry Smith   } else if (oshift == 1) {
593b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
603b2fbd54SBarry Smith     int nz = a->i[n] + 1;
613b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]++;
623b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]++;
633b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
643b2fbd54SBarry Smith   } else {
653b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
663b2fbd54SBarry Smith   }
673b2fbd54SBarry Smith 
683a40ed3dSBarry Smith   PetscFunctionReturn(0);
693b2fbd54SBarry Smith }
703b2fbd54SBarry Smith 
715615d1e5SSatish Balay #undef __FUNC__
72d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ"
733b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
743b2fbd54SBarry Smith                                 PetscTruth *done)
753b2fbd54SBarry Smith {
763b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
773b2fbd54SBarry Smith   int         i,n = a->mbs;
783b2fbd54SBarry Smith 
793a40ed3dSBarry Smith   PetscFunctionBegin;
803a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
813b2fbd54SBarry Smith   if (symmetric) {
823b2fbd54SBarry Smith     PetscFree(*ia);
833b2fbd54SBarry Smith     PetscFree(*ja);
84af5da2bfSBarry Smith   } else if (oshift == 1) {
853b2fbd54SBarry Smith     int nz = a->i[n];
863b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]--;
873b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]--;
883b2fbd54SBarry Smith   }
893a40ed3dSBarry Smith   PetscFunctionReturn(0);
903b2fbd54SBarry Smith }
913b2fbd54SBarry Smith 
922d61bbb3SSatish Balay #undef __FUNC__
932d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ"
942d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs)
952d61bbb3SSatish Balay {
962d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data;
972d61bbb3SSatish Balay 
982d61bbb3SSatish Balay   PetscFunctionBegin;
992d61bbb3SSatish Balay   *bs = baij->bs;
1002d61bbb3SSatish Balay   PetscFunctionReturn(0);
1012d61bbb3SSatish Balay }
1022d61bbb3SSatish Balay 
1032d61bbb3SSatish Balay 
1042d61bbb3SSatish Balay #undef __FUNC__
1052d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ"
106e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
1072d61bbb3SSatish Balay {
1082d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
109e51c0b9cSSatish Balay   int         ierr;
1102d61bbb3SSatish Balay 
11194d884c6SBarry Smith   if (--A->refct > 0) PetscFunctionReturn(0);
11294d884c6SBarry Smith 
11394d884c6SBarry Smith   if (A->mapping) {
11494d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr);
11594d884c6SBarry Smith   }
11694d884c6SBarry Smith   if (A->bmapping) {
11794d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr);
11894d884c6SBarry Smith   }
11961b13de0SBarry Smith   if (A->rmap) {
12061b13de0SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
12161b13de0SBarry Smith   }
12261b13de0SBarry Smith   if (A->cmap) {
12361b13de0SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
12461b13de0SBarry Smith   }
1252d61bbb3SSatish Balay #if defined(USE_PETSC_LOG)
126e1311b90SBarry Smith   PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
1272d61bbb3SSatish Balay #endif
1282d61bbb3SSatish Balay   PetscFree(a->a);
1292d61bbb3SSatish Balay   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
1302d61bbb3SSatish Balay   if (a->diag) PetscFree(a->diag);
1312d61bbb3SSatish Balay   if (a->ilen) PetscFree(a->ilen);
1322d61bbb3SSatish Balay   if (a->imax) PetscFree(a->imax);
1332d61bbb3SSatish Balay   if (a->solve_work) PetscFree(a->solve_work);
1342d61bbb3SSatish Balay   if (a->mult_work) PetscFree(a->mult_work);
135e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
1362d61bbb3SSatish Balay   PetscFree(a);
1372d61bbb3SSatish Balay   PLogObjectDestroy(A);
1382d61bbb3SSatish Balay   PetscHeaderDestroy(A);
1392d61bbb3SSatish Balay   PetscFunctionReturn(0);
1402d61bbb3SSatish Balay }
1412d61bbb3SSatish Balay 
1422d61bbb3SSatish Balay #undef __FUNC__
1432d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ"
1442d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1452d61bbb3SSatish Balay {
1462d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1472d61bbb3SSatish Balay 
1482d61bbb3SSatish Balay   PetscFunctionBegin;
1492d61bbb3SSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
1502d61bbb3SSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
1512d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
1522d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
1532d61bbb3SSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
1542d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
1552d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew       = -2;
1562d61bbb3SSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
1572d61bbb3SSatish Balay   else if (op == MAT_ROWS_SORTED ||
1582d61bbb3SSatish Balay            op == MAT_ROWS_UNSORTED ||
1592d61bbb3SSatish Balay            op == MAT_SYMMETRIC ||
1602d61bbb3SSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
1612d61bbb3SSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
162b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
163b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE) {
164981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
1652d61bbb3SSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
1662d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
1672d61bbb3SSatish Balay   } else {
1682d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
1692d61bbb3SSatish Balay   }
1702d61bbb3SSatish Balay   PetscFunctionReturn(0);
1712d61bbb3SSatish Balay }
1722d61bbb3SSatish Balay 
1732d61bbb3SSatish Balay 
1742d61bbb3SSatish Balay #undef __FUNC__
1752d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ"
1762d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
1772d61bbb3SSatish Balay {
1782d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1792d61bbb3SSatish Balay 
1802d61bbb3SSatish Balay   PetscFunctionBegin;
1812d61bbb3SSatish Balay   if (m) *m = a->m;
1822d61bbb3SSatish Balay   if (n) *n = a->n;
1832d61bbb3SSatish Balay   PetscFunctionReturn(0);
1842d61bbb3SSatish Balay }
1852d61bbb3SSatish Balay 
1862d61bbb3SSatish Balay #undef __FUNC__
1872d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ"
1882d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
1892d61bbb3SSatish Balay {
1902d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1912d61bbb3SSatish Balay 
1922d61bbb3SSatish Balay   PetscFunctionBegin;
1932d61bbb3SSatish Balay   *m = 0; *n = a->m;
1942d61bbb3SSatish Balay   PetscFunctionReturn(0);
1952d61bbb3SSatish Balay }
1962d61bbb3SSatish Balay 
1972d61bbb3SSatish Balay #undef __FUNC__
1982d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ"
1992d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2002d61bbb3SSatish Balay {
2012d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2022d61bbb3SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
2032d61bbb3SSatish Balay   Scalar      *aa,*v_i,*aa_i;
2042d61bbb3SSatish Balay 
2052d61bbb3SSatish Balay   PetscFunctionBegin;
2062d61bbb3SSatish Balay   bs  = a->bs;
2072d61bbb3SSatish Balay   ai  = a->i;
2082d61bbb3SSatish Balay   aj  = a->j;
2092d61bbb3SSatish Balay   aa  = a->a;
2102d61bbb3SSatish Balay   bs2 = a->bs2;
2112d61bbb3SSatish Balay 
2122d61bbb3SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
2132d61bbb3SSatish Balay 
2142d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2152d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2162d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2172d61bbb3SSatish Balay   *nz = bs*M;
2182d61bbb3SSatish Balay 
2192d61bbb3SSatish Balay   if (v) {
2202d61bbb3SSatish Balay     *v = 0;
2212d61bbb3SSatish Balay     if (*nz) {
2222d61bbb3SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
2232d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2242d61bbb3SSatish Balay         v_i  = *v + i*bs;
2252d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2262d61bbb3SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
2272d61bbb3SSatish Balay       }
2282d61bbb3SSatish Balay     }
2292d61bbb3SSatish Balay   }
2302d61bbb3SSatish Balay 
2312d61bbb3SSatish Balay   if (idx) {
2322d61bbb3SSatish Balay     *idx = 0;
2332d61bbb3SSatish Balay     if (*nz) {
2342d61bbb3SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
2352d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2362d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2372d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2382d61bbb3SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
2392d61bbb3SSatish Balay       }
2402d61bbb3SSatish Balay     }
2412d61bbb3SSatish Balay   }
2422d61bbb3SSatish Balay   PetscFunctionReturn(0);
2432d61bbb3SSatish Balay }
2442d61bbb3SSatish Balay 
2452d61bbb3SSatish Balay #undef __FUNC__
2462d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ"
2472d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2482d61bbb3SSatish Balay {
2492d61bbb3SSatish Balay   PetscFunctionBegin;
2502d61bbb3SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
2512d61bbb3SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
2522d61bbb3SSatish Balay   PetscFunctionReturn(0);
2532d61bbb3SSatish Balay }
2542d61bbb3SSatish Balay 
2552d61bbb3SSatish Balay #undef __FUNC__
2562d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ"
2572d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
2582d61bbb3SSatish Balay {
2592d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
2602d61bbb3SSatish Balay   Mat         C;
2612d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
2622d61bbb3SSatish Balay   int         *rows,*cols,bs2=a->bs2;
2632d61bbb3SSatish Balay   Scalar      *array=a->a;
2642d61bbb3SSatish Balay 
2652d61bbb3SSatish Balay   PetscFunctionBegin;
2662d61bbb3SSatish Balay   if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
2672d61bbb3SSatish Balay   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
2682d61bbb3SSatish Balay   PetscMemzero(col,(1+nbs)*sizeof(int));
2692d61bbb3SSatish Balay 
2702d61bbb3SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
2712d61bbb3SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
2722d61bbb3SSatish Balay   PetscFree(col);
2732d61bbb3SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
2742d61bbb3SSatish Balay   cols = rows + bs;
2752d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2762d61bbb3SSatish Balay     cols[0] = i*bs;
2772d61bbb3SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
2782d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
2792d61bbb3SSatish Balay     for ( j=0; j<len; j++ ) {
2802d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
2812d61bbb3SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
2822d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
2832d61bbb3SSatish Balay       array += bs2;
2842d61bbb3SSatish Balay     }
2852d61bbb3SSatish Balay   }
2862d61bbb3SSatish Balay   PetscFree(rows);
2872d61bbb3SSatish Balay 
2882d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2892d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
2902d61bbb3SSatish Balay 
2912d61bbb3SSatish Balay   if (B != PETSC_NULL) {
2922d61bbb3SSatish Balay     *B = C;
2932d61bbb3SSatish Balay   } else {
294f830108cSBarry Smith     PetscOps *Abops;
295cc2dc46cSBarry Smith     MatOps   Aops;
296f830108cSBarry Smith 
2972d61bbb3SSatish Balay     /* This isn't really an in-place transpose */
2982d61bbb3SSatish Balay     PetscFree(a->a);
2992d61bbb3SSatish Balay     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
3002d61bbb3SSatish Balay     if (a->diag) PetscFree(a->diag);
3012d61bbb3SSatish Balay     if (a->ilen) PetscFree(a->ilen);
3022d61bbb3SSatish Balay     if (a->imax) PetscFree(a->imax);
3032d61bbb3SSatish Balay     if (a->solve_work) PetscFree(a->solve_work);
3042d61bbb3SSatish Balay     PetscFree(a);
305f830108cSBarry Smith 
3067413bad6SBarry Smith 
3077413bad6SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
3087413bad6SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
3097413bad6SBarry Smith 
310f830108cSBarry Smith     /*
311f830108cSBarry Smith        This is horrible, horrible code. We need to keep the
312f830108cSBarry Smith       A pointers for the bops and ops but copy everything
313f830108cSBarry Smith       else from C.
314f830108cSBarry Smith     */
315f830108cSBarry Smith     Abops    = A->bops;
316f830108cSBarry Smith     Aops     = A->ops;
3172d61bbb3SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
318f830108cSBarry Smith     A->bops  = Abops;
319f830108cSBarry Smith     A->ops   = Aops;
32027a8da17SBarry Smith     A->qlist = 0;
321f830108cSBarry Smith 
3222d61bbb3SSatish Balay     PetscHeaderDestroy(C);
3232d61bbb3SSatish Balay   }
3242d61bbb3SSatish Balay   PetscFunctionReturn(0);
3252d61bbb3SSatish Balay }
3262d61bbb3SSatish Balay 
3272d61bbb3SSatish Balay 
3282d61bbb3SSatish Balay 
3293b2fbd54SBarry Smith 
3305615d1e5SSatish Balay #undef __FUNC__
331d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary"
332b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
3332593348eSBarry Smith {
334b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3359df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
336b6490206SBarry Smith   Scalar      *aa;
3372593348eSBarry Smith 
3383a40ed3dSBarry Smith   PetscFunctionBegin;
33990ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
3402593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
3412593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3423b2fbd54SBarry Smith 
3432593348eSBarry Smith   col_lens[1] = a->m;
3442593348eSBarry Smith   col_lens[2] = a->n;
3457e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3462593348eSBarry Smith 
3472593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
348b6490206SBarry Smith   count = 0;
349b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
350b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
351b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
352b6490206SBarry Smith     }
3532593348eSBarry Smith   }
3540752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3552593348eSBarry Smith   PetscFree(col_lens);
3562593348eSBarry Smith 
3572593348eSBarry Smith   /* store column indices (zero start index) */
35866963ce1SSatish Balay   jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj);
359b6490206SBarry Smith   count = 0;
360b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
361b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
362b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
363b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
364b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3652593348eSBarry Smith         }
3662593348eSBarry Smith       }
367b6490206SBarry Smith     }
368b6490206SBarry Smith   }
3690752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr);
370b6490206SBarry Smith   PetscFree(jj);
3712593348eSBarry Smith 
3722593348eSBarry Smith   /* store nonzero values */
37366963ce1SSatish Balay   aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa);
374b6490206SBarry Smith   count = 0;
375b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
376b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
377b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
378b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
3797e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
380b6490206SBarry Smith         }
381b6490206SBarry Smith       }
382b6490206SBarry Smith     }
383b6490206SBarry Smith   }
3840752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
385b6490206SBarry Smith   PetscFree(aa);
3863a40ed3dSBarry Smith   PetscFunctionReturn(0);
3872593348eSBarry Smith }
3882593348eSBarry Smith 
3895615d1e5SSatish Balay #undef __FUNC__
390d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII"
391b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
3922593348eSBarry Smith {
393b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3949df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
3952593348eSBarry Smith   FILE        *fd;
3962593348eSBarry Smith   char        *outputname;
3972593348eSBarry Smith 
3983a40ed3dSBarry Smith   PetscFunctionBegin;
39990ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
4002593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
40190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
402639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
4037ddc982cSLois Curfman McInnes     fprintf(fd,"  block size is %d\n",bs);
4042593348eSBarry Smith   }
405639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
406a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
4072593348eSBarry Smith   }
408639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
40944cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
41044cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
41144cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
41244cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
41344cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
4143a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
415e20fef11SSatish Balay           if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0)
41644cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
417e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
418e20fef11SSatish Balay           else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0)
419766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
420e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));
421e20fef11SSatish Balay           else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0)
422e20fef11SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
42344cd7ae7SLois Curfman McInnes #else
42444cd7ae7SLois Curfman McInnes           if (a->a[bs2*k + l*bs + j] != 0.0)
42544cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
42644cd7ae7SLois Curfman McInnes #endif
42744cd7ae7SLois Curfman McInnes           }
42844cd7ae7SLois Curfman McInnes         }
42944cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
43044cd7ae7SLois Curfman McInnes       }
43144cd7ae7SLois Curfman McInnes     }
43244cd7ae7SLois Curfman McInnes   }
4332593348eSBarry Smith   else {
434b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
435b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
436b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
437b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
438b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
4393a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
440e20fef11SSatish Balay           if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) {
44188685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
442e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
44388685aaeSLois Curfman McInnes           }
444e20fef11SSatish Balay           else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) {
445766eeae4SLois Curfman McInnes             fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
446e20fef11SSatish Balay               PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));
447766eeae4SLois Curfman McInnes           }
44888685aaeSLois Curfman McInnes           else {
449e20fef11SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
45088685aaeSLois Curfman McInnes           }
45188685aaeSLois Curfman McInnes #else
4527e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
45388685aaeSLois Curfman McInnes #endif
4542593348eSBarry Smith           }
4552593348eSBarry Smith         }
4562593348eSBarry Smith         fprintf(fd,"\n");
4572593348eSBarry Smith       }
4582593348eSBarry Smith     }
459b6490206SBarry Smith   }
4602593348eSBarry Smith   fflush(fd);
4613a40ed3dSBarry Smith   PetscFunctionReturn(0);
4622593348eSBarry Smith }
4632593348eSBarry Smith 
4645615d1e5SSatish Balay #undef __FUNC__
465d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw"
4663270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
4673270192aSSatish Balay {
4683270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
4693270192aSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2;
4703270192aSSatish Balay   double       xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
4713270192aSSatish Balay   Scalar       *aa;
4723270192aSSatish Balay   Draw         draw;
4733270192aSSatish Balay   DrawButton   button;
4743270192aSSatish Balay   PetscTruth   isnull;
4753270192aSSatish Balay 
4763a40ed3dSBarry Smith   PetscFunctionBegin;
4773a40ed3dSBarry Smith   ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr);
4783a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
4793270192aSSatish Balay 
4803270192aSSatish Balay   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
4813270192aSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
4823270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
4833270192aSSatish Balay   /* loop over matrix elements drawing boxes */
4843270192aSSatish Balay   color = DRAW_BLUE;
4853270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
4863270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
4873270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
4883270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
4893270192aSSatish Balay       aa = a->a + j*bs2;
4903270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
4913270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
4923270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
493b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
4943270192aSSatish Balay         }
4953270192aSSatish Balay       }
4963270192aSSatish Balay     }
4973270192aSSatish Balay   }
4983270192aSSatish Balay   color = DRAW_CYAN;
4993270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5003270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5013270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5023270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5033270192aSSatish Balay       aa = a->a + j*bs2;
5043270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5053270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5063270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
507b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5083270192aSSatish Balay         }
5093270192aSSatish Balay       }
5103270192aSSatish Balay     }
5113270192aSSatish Balay   }
5123270192aSSatish Balay 
5133270192aSSatish Balay   color = DRAW_RED;
5143270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5153270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5163270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5173270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5183270192aSSatish Balay       aa = a->a + j*bs2;
5193270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5203270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5213270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
522b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5233270192aSSatish Balay         }
5243270192aSSatish Balay       }
5253270192aSSatish Balay     }
5263270192aSSatish Balay   }
5273270192aSSatish Balay 
52855843e3eSBarry Smith   DrawSynchronizedFlush(draw);
5293270192aSSatish Balay   DrawGetPause(draw,&pause);
5303a40ed3dSBarry Smith   if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);}
5313270192aSSatish Balay 
5323270192aSSatish Balay   /* allow the matrix to zoom or shrink */
53355843e3eSBarry Smith   ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
5343270192aSSatish Balay   while (button != BUTTON_RIGHT) {
53555843e3eSBarry Smith     DrawSynchronizedClear(draw);
5363270192aSSatish Balay     if (button == BUTTON_LEFT) scale = .5;
5373270192aSSatish Balay     else if (button == BUTTON_CENTER) scale = 2.;
5383270192aSSatish Balay     xl = scale*(xl + w - xc) + xc - w*scale;
5393270192aSSatish Balay     xr = scale*(xr - w - xc) + xc + w*scale;
5403270192aSSatish Balay     yl = scale*(yl + h - yc) + yc - h*scale;
5413270192aSSatish Balay     yr = scale*(yr - h - yc) + yc + h*scale;
5423270192aSSatish Balay     w *= scale; h *= scale;
5433270192aSSatish Balay     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
5443270192aSSatish Balay     color = DRAW_BLUE;
5453270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5463270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5473270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5483270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5493270192aSSatish Balay         aa = a->a + j*bs2;
5503270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5513270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5523270192aSSatish Balay             if (PetscReal(*aa++) >=  0.) continue;
553b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5543270192aSSatish Balay           }
5553270192aSSatish Balay         }
5563270192aSSatish Balay       }
5573270192aSSatish Balay     }
5583270192aSSatish Balay     color = DRAW_CYAN;
5593270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5603270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5613270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5623270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5633270192aSSatish Balay         aa = a->a + j*bs2;
5643270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5653270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5663270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
567b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5683270192aSSatish Balay           }
5693270192aSSatish Balay         }
5703270192aSSatish Balay       }
5713270192aSSatish Balay     }
5723270192aSSatish Balay 
5733270192aSSatish Balay     color = DRAW_RED;
5743270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5753270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5763270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5773270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
5783270192aSSatish Balay         aa = a->a + j*bs2;
5793270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
5803270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
5813270192aSSatish Balay             if (PetscReal(*aa++) <= 0.) continue;
582b34f160eSSatish Balay             DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5833270192aSSatish Balay           }
5843270192aSSatish Balay         }
5853270192aSSatish Balay       }
5863270192aSSatish Balay     }
58755843e3eSBarry Smith     ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0);
5883270192aSSatish Balay   }
5893a40ed3dSBarry Smith   PetscFunctionReturn(0);
5903270192aSSatish Balay }
5913270192aSSatish Balay 
5925615d1e5SSatish Balay #undef __FUNC__
593d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ"
594e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer)
5952593348eSBarry Smith {
59619bcc07fSBarry Smith   ViewerType  vtype;
59719bcc07fSBarry Smith   int         ierr;
5982593348eSBarry Smith 
5993a40ed3dSBarry Smith   PetscFunctionBegin;
60019bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
60119bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
602a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported");
6033a40ed3dSBarry Smith   } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
6043a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
6053a40ed3dSBarry Smith   } else if (vtype == BINARY_FILE_VIEWER) {
6063a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
6073a40ed3dSBarry Smith   } else if (vtype == DRAW_VIEWER) {
6083a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
6095cd90555SBarry Smith   } else {
6105cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
6112593348eSBarry Smith   }
6123a40ed3dSBarry Smith   PetscFunctionReturn(0);
6132593348eSBarry Smith }
614b6490206SBarry Smith 
615cd0e1443SSatish Balay 
6165615d1e5SSatish Balay #undef __FUNC__
6172d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ"
6182d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
619cd0e1443SSatish Balay {
620cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6212d61bbb3SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
6222d61bbb3SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
6232d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6242d61bbb3SSatish Balay   Scalar     *ap, *aa = a->a, zero = 0.0;
625cd0e1443SSatish Balay 
6263a40ed3dSBarry Smith   PetscFunctionBegin;
6272d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
628cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
629a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
630a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
6312d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6322c3acbe9SBarry Smith     nrow = ailen[brow];
6332d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
634a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
635a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6362d61bbb3SSatish Balay       col  = in[l] ;
6372d61bbb3SSatish Balay       bcol = col/bs;
6382d61bbb3SSatish Balay       cidx = col%bs;
6392d61bbb3SSatish Balay       ridx = row%bs;
6402d61bbb3SSatish Balay       high = nrow;
6412d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6422d61bbb3SSatish Balay       while (high-low > 5) {
643cd0e1443SSatish Balay         t = (low+high)/2;
644cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
645cd0e1443SSatish Balay         else             low  = t;
646cd0e1443SSatish Balay       }
647cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
648cd0e1443SSatish Balay         if (rp[i] > bcol) break;
649cd0e1443SSatish Balay         if (rp[i] == bcol) {
6502d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6512d61bbb3SSatish Balay           goto finished;
652cd0e1443SSatish Balay         }
653cd0e1443SSatish Balay       }
6542d61bbb3SSatish Balay       *v++ = zero;
6552d61bbb3SSatish Balay       finished:;
656cd0e1443SSatish Balay     }
657cd0e1443SSatish Balay   }
6583a40ed3dSBarry Smith   PetscFunctionReturn(0);
659cd0e1443SSatish Balay }
660cd0e1443SSatish Balay 
6612d61bbb3SSatish Balay 
6625615d1e5SSatish Balay #undef __FUNC__
66305a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ"
66492c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
66592c4ed94SBarry Smith {
66692c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6678a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
668dd9472c6SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented;
669dd9472c6SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
6700e324ae4SSatish Balay   Scalar      *ap,*value=v,*aa=a->a,*bap;
67192c4ed94SBarry Smith 
6723a40ed3dSBarry Smith   PetscFunctionBegin;
6730e324ae4SSatish Balay   if (roworiented) {
6740e324ae4SSatish Balay     stepval = (n-1)*bs;
6750e324ae4SSatish Balay   } else {
6760e324ae4SSatish Balay     stepval = (m-1)*bs;
6770e324ae4SSatish Balay   }
67892c4ed94SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
67992c4ed94SBarry Smith     row  = im[k];
6803a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
681a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
682a8c6a408SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
68392c4ed94SBarry Smith #endif
68492c4ed94SBarry Smith     rp   = aj + ai[row];
68592c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
68692c4ed94SBarry Smith     rmax = imax[row];
68792c4ed94SBarry Smith     nrow = ailen[row];
68892c4ed94SBarry Smith     low  = 0;
68992c4ed94SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
6903a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
691a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
692a8c6a408SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
69392c4ed94SBarry Smith #endif
69492c4ed94SBarry Smith       col = in[l];
69592c4ed94SBarry Smith       if (roworiented) {
6960e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
6970e324ae4SSatish Balay       } else {
6980e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
69992c4ed94SBarry Smith       }
70092c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
70192c4ed94SBarry Smith       while (high-low > 7) {
70292c4ed94SBarry Smith         t = (low+high)/2;
70392c4ed94SBarry Smith         if (rp[t] > col) high = t;
70492c4ed94SBarry Smith         else             low  = t;
70592c4ed94SBarry Smith       }
70692c4ed94SBarry Smith       for ( i=low; i<high; i++ ) {
70792c4ed94SBarry Smith         if (rp[i] > col) break;
70892c4ed94SBarry Smith         if (rp[i] == col) {
7098a84c255SSatish Balay           bap  = ap +  bs2*i;
7100e324ae4SSatish Balay           if (roworiented) {
7118a84c255SSatish Balay             if (is == ADD_VALUES) {
712dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
713dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7148a84c255SSatish Balay                   bap[jj] += *value++;
715dd9472c6SBarry Smith                 }
716dd9472c6SBarry Smith               }
7170e324ae4SSatish Balay             } else {
718dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
719dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7200e324ae4SSatish Balay                   bap[jj] = *value++;
7218a84c255SSatish Balay                 }
722dd9472c6SBarry Smith               }
723dd9472c6SBarry Smith             }
7240e324ae4SSatish Balay           } else {
7250e324ae4SSatish Balay             if (is == ADD_VALUES) {
726dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
727dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7280e324ae4SSatish Balay                   *bap++ += *value++;
729dd9472c6SBarry Smith                 }
730dd9472c6SBarry Smith               }
7310e324ae4SSatish Balay             } else {
732dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
733dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7340e324ae4SSatish Balay                   *bap++  = *value++;
7350e324ae4SSatish Balay                 }
7368a84c255SSatish Balay               }
737dd9472c6SBarry Smith             }
738dd9472c6SBarry Smith           }
739f1241b54SBarry Smith           goto noinsert2;
74092c4ed94SBarry Smith         }
74192c4ed94SBarry Smith       }
74289280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
743a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
74492c4ed94SBarry Smith       if (nrow >= rmax) {
74592c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
74692c4ed94SBarry Smith         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
74792c4ed94SBarry Smith         Scalar *new_a;
74892c4ed94SBarry Smith 
749a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
75089280ab3SLois Curfman McInnes 
75192c4ed94SBarry Smith         /* malloc new storage space */
75292c4ed94SBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
75392c4ed94SBarry Smith         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
75492c4ed94SBarry Smith         new_j   = (int *) (new_a + bs2*new_nz);
75592c4ed94SBarry Smith         new_i   = new_j + new_nz;
75692c4ed94SBarry Smith 
75792c4ed94SBarry Smith         /* copy over old data into new slots */
75892c4ed94SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
75992c4ed94SBarry Smith         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
76092c4ed94SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));
76192c4ed94SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow);
762dd9472c6SBarry Smith         PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));
76392c4ed94SBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(Scalar));
76492c4ed94SBarry Smith         PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
765dd9472c6SBarry Smith         PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(Scalar));
76692c4ed94SBarry Smith         /* free up old matrix storage */
76792c4ed94SBarry Smith         PetscFree(a->a);
76892c4ed94SBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
76992c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
77092c4ed94SBarry Smith         a->singlemalloc = 1;
77192c4ed94SBarry Smith 
77292c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
77392c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
77492c4ed94SBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
77592c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
77692c4ed94SBarry Smith         a->reallocs++;
77792c4ed94SBarry Smith         a->nz++;
77892c4ed94SBarry Smith       }
77992c4ed94SBarry Smith       N = nrow++ - 1;
78092c4ed94SBarry Smith       /* shift up all the later entries in this row */
78192c4ed94SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
78292c4ed94SBarry Smith         rp[ii+1] = rp[ii];
78392c4ed94SBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
78492c4ed94SBarry Smith       }
78592c4ed94SBarry Smith       if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
78692c4ed94SBarry Smith       rp[i] = col;
7878a84c255SSatish Balay       bap   = ap +  bs2*i;
7880e324ae4SSatish Balay       if (roworiented) {
789dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval) {
790dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs ) {
7910e324ae4SSatish Balay             bap[jj] = *value++;
792dd9472c6SBarry Smith           }
793dd9472c6SBarry Smith         }
7940e324ae4SSatish Balay       } else {
795dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval ) {
796dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++ ) {
7970e324ae4SSatish Balay             *bap++  = *value++;
7980e324ae4SSatish Balay           }
799dd9472c6SBarry Smith         }
800dd9472c6SBarry Smith       }
801f1241b54SBarry Smith       noinsert2:;
80292c4ed94SBarry Smith       low = i;
80392c4ed94SBarry Smith     }
80492c4ed94SBarry Smith     ailen[row] = nrow;
80592c4ed94SBarry Smith   }
8063a40ed3dSBarry Smith   PetscFunctionReturn(0);
80792c4ed94SBarry Smith }
80892c4ed94SBarry Smith 
809f2501298SSatish Balay 
8105615d1e5SSatish Balay #undef __FUNC__
8115615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ"
8128f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
813584200bdSSatish Balay {
814584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
815584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
816a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
81743ee02c3SBarry Smith   int        mbs = a->mbs, bs2 = a->bs2,rmax = 0;
818584200bdSSatish Balay   Scalar     *aa = a->a, *ap;
819584200bdSSatish Balay 
8203a40ed3dSBarry Smith   PetscFunctionBegin;
8213a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
822584200bdSSatish Balay 
82343ee02c3SBarry Smith   if (m) rmax = ailen[0];
824584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
825584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
826584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
827d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
828584200bdSSatish Balay     if (fshift) {
829a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
830584200bdSSatish Balay       N = ailen[i];
831584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
832584200bdSSatish Balay         ip[j-fshift] = ip[j];
8337e67e3f9SSatish Balay         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar));
834584200bdSSatish Balay       }
835584200bdSSatish Balay     }
836584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
837584200bdSSatish Balay   }
838584200bdSSatish Balay   if (mbs) {
839584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
840584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
841584200bdSSatish Balay   }
842584200bdSSatish Balay   /* reset ilen and imax for each row */
843584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
844584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
845584200bdSSatish Balay   }
846a7c10996SSatish Balay   a->nz = ai[mbs];
847584200bdSSatish Balay 
848584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
849584200bdSSatish Balay   if (fshift && a->diag) {
850584200bdSSatish Balay     PetscFree(a->diag);
851584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
852584200bdSSatish Balay     a->diag = 0;
853584200bdSSatish Balay   }
8543d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
855219d9a1aSLois Curfman McInnes            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
8563d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
857584200bdSSatish Balay            a->reallocs);
858d402145bSBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
859e2f3b5e9SSatish Balay   a->reallocs          = 0;
8604e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
8614e220ebcSLois Curfman McInnes 
8623a40ed3dSBarry Smith   PetscFunctionReturn(0);
863584200bdSSatish Balay }
864584200bdSSatish Balay 
8655a838052SSatish Balay 
866d9b7c43dSSatish Balay /* idx should be of length atlease bs */
8675615d1e5SSatish Balay #undef __FUNC__
8685615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block"
869d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg)
870d9b7c43dSSatish Balay {
871d9b7c43dSSatish Balay   int i,row;
8723a40ed3dSBarry Smith 
8733a40ed3dSBarry Smith   PetscFunctionBegin;
874d9b7c43dSSatish Balay   row = idx[0];
8753a40ed3dSBarry Smith   if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
876d9b7c43dSSatish Balay 
877d9b7c43dSSatish Balay   for ( i=1; i<bs; i++ ) {
8783a40ed3dSBarry Smith     if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); }
879d9b7c43dSSatish Balay   }
880d9b7c43dSSatish Balay   *flg = PETSC_TRUE;
8813a40ed3dSBarry Smith   PetscFunctionReturn(0);
882d9b7c43dSSatish Balay }
883d9b7c43dSSatish Balay 
8845615d1e5SSatish Balay #undef __FUNC__
8855615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ"
8868f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
887d9b7c43dSSatish Balay {
888d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
889d9b7c43dSSatish Balay   IS          is_local;
890d9b7c43dSSatish Balay   int         ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2;
891d9b7c43dSSatish Balay   PetscTruth  flg;
892d9b7c43dSSatish Balay   Scalar      zero = 0.0,*aa;
893d9b7c43dSSatish Balay 
8943a40ed3dSBarry Smith   PetscFunctionBegin;
895d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
896d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
897d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
898537820f0SBarry Smith   ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr);
899d9b7c43dSSatish Balay   ierr = ISSort(is_local); CHKERRQ(ierr);
900d9b7c43dSSatish Balay   ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr);
901d9b7c43dSSatish Balay 
902d9b7c43dSSatish Balay   i = 0;
903d9b7c43dSSatish Balay   while (i < is_n) {
904a8c6a408SBarry Smith     if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range");
905d9b7c43dSSatish Balay     flg = PETSC_FALSE;
906d9b7c43dSSatish Balay     if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); }
907d9b7c43dSSatish Balay     count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs;
908d9b7c43dSSatish Balay     aa    = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs);
9092d61bbb3SSatish Balay     if (flg) { /* There exists a block of rows to be Zerowed */
910a07cd24cSSatish Balay       if (baij->ilen[rows[i]/bs] > 0) {
9112d61bbb3SSatish Balay         PetscMemzero(aa,count*bs*sizeof(Scalar));
9122d61bbb3SSatish Balay         baij->ilen[rows[i]/bs] = 1;
9132d61bbb3SSatish Balay         baij->j[baij->i[rows[i]/bs]] = rows[i]/bs;
914a07cd24cSSatish Balay       }
9152d61bbb3SSatish Balay       i += bs;
9162d61bbb3SSatish Balay     } else { /* Zero out only the requested row */
917d9b7c43dSSatish Balay       for ( j=0; j<count; j++ ) {
918d9b7c43dSSatish Balay         aa[0] = zero;
919d9b7c43dSSatish Balay         aa+=bs;
920d9b7c43dSSatish Balay       }
921d9b7c43dSSatish Balay       i++;
922d9b7c43dSSatish Balay     }
923d9b7c43dSSatish Balay   }
924d9b7c43dSSatish Balay   if (diag) {
925d9b7c43dSSatish Balay     for ( j=0; j<is_n; j++ ) {
926f830108cSBarry Smith       ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
9272d61bbb3SSatish Balay       /* ierr = MatSetValues(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); */
928d9b7c43dSSatish Balay     }
929d9b7c43dSSatish Balay   }
930d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr);
931d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr);
932d9b7c43dSSatish Balay   ierr = ISDestroy(is_local); CHKERRQ(ierr);
9339a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9343a40ed3dSBarry Smith   PetscFunctionReturn(0);
935d9b7c43dSSatish Balay }
9361c351548SSatish Balay 
9375615d1e5SSatish Balay #undef __FUNC__
9382d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ"
9392d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
9402d61bbb3SSatish Balay {
9412d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9422d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
9432d61bbb3SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
9442d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
9452d61bbb3SSatish Balay   int          ridx,cidx,bs2=a->bs2;
9462d61bbb3SSatish Balay   Scalar      *ap,value,*aa=a->a,*bap;
9472d61bbb3SSatish Balay 
9482d61bbb3SSatish Balay   PetscFunctionBegin;
9492d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
9502d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
9512d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9522d61bbb3SSatish Balay     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
9532d61bbb3SSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
9542d61bbb3SSatish Balay #endif
9552d61bbb3SSatish Balay     rp   = aj + ai[brow];
9562d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
9572d61bbb3SSatish Balay     rmax = imax[brow];
9582d61bbb3SSatish Balay     nrow = ailen[brow];
9592d61bbb3SSatish Balay     low  = 0;
9602d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
9612d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9622d61bbb3SSatish Balay       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
9632d61bbb3SSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
9642d61bbb3SSatish Balay #endif
9652d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
9662d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
9672d61bbb3SSatish Balay       if (roworiented) {
9682d61bbb3SSatish Balay         value = *v++;
9692d61bbb3SSatish Balay       } else {
9702d61bbb3SSatish Balay         value = v[k + l*m];
9712d61bbb3SSatish Balay       }
9722d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
9732d61bbb3SSatish Balay       while (high-low > 7) {
9742d61bbb3SSatish Balay         t = (low+high)/2;
9752d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
9762d61bbb3SSatish Balay         else              low  = t;
9772d61bbb3SSatish Balay       }
9782d61bbb3SSatish Balay       for ( i=low; i<high; i++ ) {
9792d61bbb3SSatish Balay         if (rp[i] > bcol) break;
9802d61bbb3SSatish Balay         if (rp[i] == bcol) {
9812d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
9822d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
9832d61bbb3SSatish Balay           else                  *bap  = value;
9842d61bbb3SSatish Balay           goto noinsert1;
9852d61bbb3SSatish Balay         }
9862d61bbb3SSatish Balay       }
9872d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
9882d61bbb3SSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9892d61bbb3SSatish Balay       if (nrow >= rmax) {
9902d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
9912d61bbb3SSatish Balay         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
9922d61bbb3SSatish Balay         Scalar *new_a;
9932d61bbb3SSatish Balay 
9942d61bbb3SSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
9952d61bbb3SSatish Balay 
9962d61bbb3SSatish Balay         /* Malloc new storage space */
9972d61bbb3SSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
9982d61bbb3SSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
9992d61bbb3SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
10002d61bbb3SSatish Balay         new_i   = new_j + new_nz;
10012d61bbb3SSatish Balay 
10022d61bbb3SSatish Balay         /* copy over old data into new slots */
10032d61bbb3SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
10042d61bbb3SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
10052d61bbb3SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
10062d61bbb3SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
10072d61bbb3SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
10082d61bbb3SSatish Balay                                                            len*sizeof(int));
10092d61bbb3SSatish Balay         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar));
10102d61bbb3SSatish Balay         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
10112d61bbb3SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
10122d61bbb3SSatish Balay                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar));
10132d61bbb3SSatish Balay         /* free up old matrix storage */
10142d61bbb3SSatish Balay         PetscFree(a->a);
10152d61bbb3SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
10162d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
10172d61bbb3SSatish Balay         a->singlemalloc = 1;
10182d61bbb3SSatish Balay 
10192d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10202d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
10212d61bbb3SSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
10222d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10232d61bbb3SSatish Balay         a->reallocs++;
10242d61bbb3SSatish Balay         a->nz++;
10252d61bbb3SSatish Balay       }
10262d61bbb3SSatish Balay       N = nrow++ - 1;
10272d61bbb3SSatish Balay       /* shift up all the later entries in this row */
10282d61bbb3SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
10292d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
10302d61bbb3SSatish Balay         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
10312d61bbb3SSatish Balay       }
10322d61bbb3SSatish Balay       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
10332d61bbb3SSatish Balay       rp[i]                      = bcol;
10342d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
10352d61bbb3SSatish Balay       noinsert1:;
10362d61bbb3SSatish Balay       low = i;
10372d61bbb3SSatish Balay     }
10382d61bbb3SSatish Balay     ailen[brow] = nrow;
10392d61bbb3SSatish Balay   }
10402d61bbb3SSatish Balay   PetscFunctionReturn(0);
10412d61bbb3SSatish Balay }
10422d61bbb3SSatish Balay 
10432d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
10442d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
10452d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
1046d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*);
10472d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**);
10482d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec);
10492d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
10502d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat);
10512d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *);
10522d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*);
10532d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec);
10542d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
10552d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
10562d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat);
10572d61bbb3SSatish Balay 
10582d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
10592d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
10602d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
10612d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
10622d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
10632d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
10642d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
10652d61bbb3SSatish Balay 
10662d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
10672d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
10682d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
10692d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
10702d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
10712d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
10722d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
10732d61bbb3SSatish Balay 
10742d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
10752d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
10762d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
10772d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
10782d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
10792d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
10802d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
10812d61bbb3SSatish Balay 
10822d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
10832d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
10842d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
10852d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
10862d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
10872d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
10882d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
10892d61bbb3SSatish Balay 
10902d61bbb3SSatish Balay /*
10912d61bbb3SSatish Balay      note: This can only work for identity for row and col. It would
10922d61bbb3SSatish Balay    be good to check this and otherwise generate an error.
10932d61bbb3SSatish Balay */
10942d61bbb3SSatish Balay #undef __FUNC__
10952d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ"
10962d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
10972d61bbb3SSatish Balay {
10982d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
10992d61bbb3SSatish Balay   Mat         outA;
11002d61bbb3SSatish Balay   int         ierr;
11012d61bbb3SSatish Balay 
11022d61bbb3SSatish Balay   PetscFunctionBegin;
11032d61bbb3SSatish Balay   if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
11042d61bbb3SSatish Balay 
11052d61bbb3SSatish Balay   outA          = inA;
11062d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
11072d61bbb3SSatish Balay   a->row        = row;
11082d61bbb3SSatish Balay   a->col        = col;
11092d61bbb3SSatish Balay 
1110e51c0b9cSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1111e51c0b9cSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
11121e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1113e51c0b9cSSatish Balay 
11142d61bbb3SSatish Balay   if (!a->solve_work) {
11152d61bbb3SSatish Balay     a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
11162d61bbb3SSatish Balay     PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
11172d61bbb3SSatish Balay   }
11182d61bbb3SSatish Balay 
11192d61bbb3SSatish Balay   if (!a->diag) {
11202d61bbb3SSatish Balay     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
11212d61bbb3SSatish Balay   }
11222d61bbb3SSatish Balay   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
11232d61bbb3SSatish Balay 
11242d61bbb3SSatish Balay   /*
11252d61bbb3SSatish Balay       Blocksize 4 has a special faster solver for ILU(0) factorization
11262d61bbb3SSatish Balay     with natural ordering
11272d61bbb3SSatish Balay   */
11282d61bbb3SSatish Balay   if (a->bs == 4) {
1129f830108cSBarry Smith     inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering;
11302d61bbb3SSatish Balay   }
11312d61bbb3SSatish Balay 
11322d61bbb3SSatish Balay   PetscFunctionReturn(0);
11332d61bbb3SSatish Balay }
11342d61bbb3SSatish Balay #undef __FUNC__
1135d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ"
1136ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1137ba4ca20aSSatish Balay {
1138ba4ca20aSSatish Balay   static int called = 0;
1139ba4ca20aSSatish Balay   MPI_Comm   comm = A->comm;
1140ba4ca20aSSatish Balay 
11413a40ed3dSBarry Smith   PetscFunctionBegin;
11423a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
114376be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");
114476be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");
11453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1146ba4ca20aSSatish Balay }
1147d9b7c43dSSatish Balay 
114827a8da17SBarry Smith #undef __FUNC__
114927a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
115027a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
115127a8da17SBarry Smith {
115227a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
115327a8da17SBarry Smith   int         i,nz,n;
115427a8da17SBarry Smith 
115527a8da17SBarry Smith   PetscFunctionBegin;
115627a8da17SBarry Smith   nz = baij->maxnz;
115727a8da17SBarry Smith   n  = baij->n;
115827a8da17SBarry Smith   for (i=0; i<nz; i++) {
115927a8da17SBarry Smith     baij->j[i] = indices[i];
116027a8da17SBarry Smith   }
116127a8da17SBarry Smith   baij->nz = nz;
116227a8da17SBarry Smith   for ( i=0; i<n; i++ ) {
116327a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
116427a8da17SBarry Smith   }
116527a8da17SBarry Smith 
116627a8da17SBarry Smith   PetscFunctionReturn(0);
116727a8da17SBarry Smith }
116827a8da17SBarry Smith 
116927a8da17SBarry Smith #undef __FUNC__
117027a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices"
117127a8da17SBarry Smith /*@
117227a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
117327a8da17SBarry Smith        in the matrix.
117427a8da17SBarry Smith 
117527a8da17SBarry Smith   Input Parameters:
117627a8da17SBarry Smith +  mat - the SeqBAIJ matrix
117727a8da17SBarry Smith -  indices - the column indices
117827a8da17SBarry Smith 
117927a8da17SBarry Smith   Notes:
118027a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
118127a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
118227a8da17SBarry Smith   of the MatSetValues() operation.
118327a8da17SBarry Smith 
118427a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
118527a8da17SBarry Smith   MatCreateSeqBAIJ().
118627a8da17SBarry Smith 
118727a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
118827a8da17SBarry Smith 
118927a8da17SBarry Smith @*/
119027a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
119127a8da17SBarry Smith {
119227a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
119327a8da17SBarry Smith 
119427a8da17SBarry Smith   PetscFunctionBegin;
119527a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
119627a8da17SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);
119727a8da17SBarry Smith          CHKERRQ(ierr);
119827a8da17SBarry Smith   if (f) {
119927a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
120027a8da17SBarry Smith   } else {
120127a8da17SBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
120227a8da17SBarry Smith   }
120327a8da17SBarry Smith   PetscFunctionReturn(0);
120427a8da17SBarry Smith }
120527a8da17SBarry Smith 
12062593348eSBarry Smith /* -------------------------------------------------------------------*/
1207cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1208cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1209cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1210cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1211cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
1212cc2dc46cSBarry Smith        MatMultTrans_SeqBAIJ,
1213cc2dc46cSBarry Smith        MatMultTransAdd_SeqBAIJ,
1214cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1215cc2dc46cSBarry Smith        0,
1216cc2dc46cSBarry Smith        0,
1217cc2dc46cSBarry Smith        0,
1218cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1219cc2dc46cSBarry Smith        0,
1220b6490206SBarry Smith        0,
1221f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1222cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1223cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1224cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1225cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1226cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1227b6490206SBarry Smith        0,
1228cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1229cc2dc46cSBarry Smith        0,
1230cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1231cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1232cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1233cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1234cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1235cc2dc46cSBarry Smith        0,
1236cc2dc46cSBarry Smith        0,
1237cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1238cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1239cc2dc46cSBarry Smith        MatGetOwnershipRange_SeqBAIJ,
1240cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1241cc2dc46cSBarry Smith        0,
1242cc2dc46cSBarry Smith        0,
1243cc2dc46cSBarry Smith        0,
1244cc2dc46cSBarry Smith        MatConvertSameType_SeqBAIJ,
1245cc2dc46cSBarry Smith        0,
1246cc2dc46cSBarry Smith        0,
1247cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1248cc2dc46cSBarry Smith        0,
1249cc2dc46cSBarry Smith        0,
1250cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1251cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1252cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1253cc2dc46cSBarry Smith        0,
1254cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1255cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1256cc2dc46cSBarry Smith        0,
1257cc2dc46cSBarry Smith        0,
1258cc2dc46cSBarry Smith        0,
1259cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
12603b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
126192c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1262cc2dc46cSBarry Smith        0,
1263cc2dc46cSBarry Smith        0,
1264cc2dc46cSBarry Smith        0,
1265cc2dc46cSBarry Smith        0,
1266cc2dc46cSBarry Smith        0,
1267cc2dc46cSBarry Smith        0,
1268d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1269cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1270cc2dc46cSBarry Smith        0,
1271cc2dc46cSBarry Smith        0,
1272cc2dc46cSBarry Smith        MatGetMaps_Petsc};
12732593348eSBarry Smith 
12745615d1e5SSatish Balay #undef __FUNC__
12755615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ"
12762593348eSBarry Smith /*@C
127744cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
127844cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
127944cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
12802bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
12812bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
12822593348eSBarry Smith 
1283db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1284db81eaa0SLois Curfman McInnes 
12852593348eSBarry Smith    Input Parameters:
1286db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
1287b6490206SBarry Smith .  bs - size of block
12882593348eSBarry Smith .  m - number of rows
12892593348eSBarry Smith .  n - number of columns
1290b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1291db81eaa0SLois Curfman McInnes -  nzz - array containing the number of block nonzeros in the various block rows
12922bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
12932593348eSBarry Smith 
12942593348eSBarry Smith    Output Parameter:
12952593348eSBarry Smith .  A - the matrix
12962593348eSBarry Smith 
12970513a670SBarry Smith    Options Database Keys:
1298db81eaa0SLois Curfman McInnes .   -mat_no_unroll - uses code that does not unroll the loops in the
1299db81eaa0SLois Curfman McInnes                      block calculations (much slower)
1300db81eaa0SLois Curfman McInnes .    -mat_block_size - size of the blocks to use
13010513a670SBarry Smith 
13022593348eSBarry Smith    Notes:
130344cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
13042593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
130544cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
13062593348eSBarry Smith 
13072593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
13082593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
13092593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
13106da5968aSLois Curfman McInnes    matrices.
13112593348eSBarry Smith 
1312db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
13132593348eSBarry Smith @*/
1314b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
13152593348eSBarry Smith {
13162593348eSBarry Smith   Mat         B;
1317b6490206SBarry Smith   Mat_SeqBAIJ *b;
13183b2fbd54SBarry Smith   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size;
13193b2fbd54SBarry Smith 
13203a40ed3dSBarry Smith   PetscFunctionBegin;
13213b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1322a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1323b6490206SBarry Smith 
13240513a670SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
13250513a670SBarry Smith 
13263a40ed3dSBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1327a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
13283a40ed3dSBarry Smith   }
13292593348eSBarry Smith 
13302593348eSBarry Smith   *A = 0;
1331f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView);
13322593348eSBarry Smith   PLogObjectCreate(B);
1333b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
133444cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
1335cc2dc46cSBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
13361a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
13371a6a6d98SBarry Smith   if (!flg) {
13387fc0212eSBarry Smith     switch (bs) {
13397fc0212eSBarry Smith     case 1:
1340f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1341f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1342f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1343f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
13447fc0212eSBarry Smith       break;
13454eeb42bcSBarry Smith     case 2:
1346f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1347f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1348f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1349f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
13506d84be18SBarry Smith       break;
1351f327dd97SBarry Smith     case 3:
1352f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1353f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1354f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1355f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
13564eeb42bcSBarry Smith       break;
13576d84be18SBarry Smith     case 4:
1358f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1359f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1360f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1361f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
13626d84be18SBarry Smith       break;
13636d84be18SBarry Smith     case 5:
1364f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1365f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1366f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1367f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
13686d84be18SBarry Smith       break;
136948e9ddb2SSatish Balay     case 7:
1370f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1371f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1372f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
137348e9ddb2SSatish Balay       break;
13747fc0212eSBarry Smith     }
13751a6a6d98SBarry Smith   }
1376e1311b90SBarry Smith   B->ops->destroy     = MatDestroy_SeqBAIJ;
1377e1311b90SBarry Smith   B->ops->view        = MatView_SeqBAIJ;
13782593348eSBarry Smith   B->factor           = 0;
13792593348eSBarry Smith   B->lupivotthreshold = 1.0;
138090f02eecSBarry Smith   B->mapping          = 0;
13812593348eSBarry Smith   b->row              = 0;
13822593348eSBarry Smith   b->col              = 0;
1383e51c0b9cSSatish Balay   b->icol             = 0;
13842593348eSBarry Smith   b->reallocs         = 0;
13852593348eSBarry Smith 
138644cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
138744cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1388a5ae1ecdSBarry Smith 
13897413bad6SBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
13907413bad6SBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1391a5ae1ecdSBarry Smith 
1392b6490206SBarry Smith   b->mbs     = mbs;
1393f2501298SSatish Balay   b->nbs     = nbs;
1394b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
13952593348eSBarry Smith   if (nnz == PETSC_NULL) {
1396119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
13972593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1398b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1399b6490206SBarry Smith     nz = nz*mbs;
14003a40ed3dSBarry Smith   } else {
14012593348eSBarry Smith     nz = 0;
1402b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
14032593348eSBarry Smith   }
14042593348eSBarry Smith 
14052593348eSBarry Smith   /* allocate the matrix space */
14067e67e3f9SSatish Balay   len   = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int);
14072593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
14087e67e3f9SSatish Balay   PetscMemzero(b->a,nz*bs2*sizeof(Scalar));
14097e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
14102593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
14112593348eSBarry Smith   b->i  = b->j + nz;
14122593348eSBarry Smith   b->singlemalloc = 1;
14132593348eSBarry Smith 
1414de6a44a3SBarry Smith   b->i[0] = 0;
1415b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
14162593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
14172593348eSBarry Smith   }
14182593348eSBarry Smith 
1419b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1420b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1421f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1422b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
14232593348eSBarry Smith 
1424b6490206SBarry Smith   b->bs               = bs;
14259df24120SSatish Balay   b->bs2              = bs2;
1426b6490206SBarry Smith   b->mbs              = mbs;
14272593348eSBarry Smith   b->nz               = 0;
142818e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
14292593348eSBarry Smith   b->sorted           = 0;
14302593348eSBarry Smith   b->roworiented      = 1;
14312593348eSBarry Smith   b->nonew            = 0;
14322593348eSBarry Smith   b->diag             = 0;
14332593348eSBarry Smith   b->solve_work       = 0;
1434de6a44a3SBarry Smith   b->mult_work        = 0;
14352593348eSBarry Smith   b->spptr            = 0;
14364e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
14374e220ebcSLois Curfman McInnes 
14382593348eSBarry Smith   *A = B;
14392593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
14402593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
144127a8da17SBarry Smith 
144227a8da17SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
144327a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
144427a8da17SBarry Smith                                      (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
144527a8da17SBarry Smith 
14463a40ed3dSBarry Smith   PetscFunctionReturn(0);
14472593348eSBarry Smith }
14482593348eSBarry Smith 
14495615d1e5SSatish Balay #undef __FUNC__
14505615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ"
1451b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
14522593348eSBarry Smith {
14532593348eSBarry Smith   Mat         C;
1454b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
145527a8da17SBarry Smith   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1456de6a44a3SBarry Smith 
14573a40ed3dSBarry Smith   PetscFunctionBegin;
1458a8c6a408SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
14592593348eSBarry Smith 
14602593348eSBarry Smith   *B = 0;
1461f830108cSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView);
14622593348eSBarry Smith   PLogObjectCreate(C);
1463b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
1464c3c66ee5SSatish Balay   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1465e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqBAIJ;
1466e1311b90SBarry Smith   C->ops->view       = MatView_SeqBAIJ;
14672593348eSBarry Smith   C->factor     = A->factor;
14682593348eSBarry Smith   c->row        = 0;
14692593348eSBarry Smith   c->col        = 0;
14702593348eSBarry Smith   C->assembled  = PETSC_TRUE;
14712593348eSBarry Smith 
147244cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
147344cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
147444cd7ae7SLois Curfman McInnes   C->M          = a->m;
147544cd7ae7SLois Curfman McInnes   C->N          = a->n;
147644cd7ae7SLois Curfman McInnes 
1477b6490206SBarry Smith   c->bs         = a->bs;
14789df24120SSatish Balay   c->bs2        = a->bs2;
1479b6490206SBarry Smith   c->mbs        = a->mbs;
1480b6490206SBarry Smith   c->nbs        = a->nbs;
14812593348eSBarry Smith 
1482b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1483b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1484b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
14852593348eSBarry Smith     c->imax[i] = a->imax[i];
14862593348eSBarry Smith     c->ilen[i] = a->ilen[i];
14872593348eSBarry Smith   }
14882593348eSBarry Smith 
14892593348eSBarry Smith   /* allocate the matrix space */
14902593348eSBarry Smith   c->singlemalloc = 1;
14917e67e3f9SSatish Balay   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int));
14922593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
14937e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1494de6a44a3SBarry Smith   c->i  = c->j + nz;
1495b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1496b6490206SBarry Smith   if (mbs > 0) {
1497de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
14982593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
14997e67e3f9SSatish Balay       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar));
15002593348eSBarry Smith     }
15012593348eSBarry Smith   }
15022593348eSBarry Smith 
1503f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
15042593348eSBarry Smith   c->sorted      = a->sorted;
15052593348eSBarry Smith   c->roworiented = a->roworiented;
15062593348eSBarry Smith   c->nonew       = a->nonew;
15072593348eSBarry Smith 
15082593348eSBarry Smith   if (a->diag) {
1509b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1510b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1511b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
15122593348eSBarry Smith       c->diag[i] = a->diag[i];
15132593348eSBarry Smith     }
15142593348eSBarry Smith   }
15152593348eSBarry Smith   else c->diag          = 0;
15162593348eSBarry Smith   c->nz                 = a->nz;
15172593348eSBarry Smith   c->maxnz              = a->maxnz;
15182593348eSBarry Smith   c->solve_work         = 0;
15192593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
15207fc0212eSBarry Smith   c->mult_work          = 0;
15212593348eSBarry Smith   *B = C;
152227a8da17SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqBAIJSetColumnIndices_C",
152327a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
152427a8da17SBarry Smith                                      (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
15253a40ed3dSBarry Smith   PetscFunctionReturn(0);
15262593348eSBarry Smith }
15272593348eSBarry Smith 
15285615d1e5SSatish Balay #undef __FUNC__
15295615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ"
153019bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
15312593348eSBarry Smith {
1532b6490206SBarry Smith   Mat_SeqBAIJ  *a;
15332593348eSBarry Smith   Mat          B;
1534de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1535b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
153635aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1537a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1538b6490206SBarry Smith   Scalar       *aa;
153919bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
15402593348eSBarry Smith 
15413a40ed3dSBarry Smith   PetscFunctionBegin;
15421a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1543de6a44a3SBarry Smith   bs2  = bs*bs;
1544b6490206SBarry Smith 
15452593348eSBarry Smith   MPI_Comm_size(comm,&size);
1546cc0fae11SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
154790ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
15480752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
1549a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
15502593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
15512593348eSBarry Smith 
1552d64ed03dSBarry Smith   if (header[3] < 0) {
1553a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1554d64ed03dSBarry Smith   }
1555d64ed03dSBarry Smith 
1556a8c6a408SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
155735aab85fSBarry Smith 
155835aab85fSBarry Smith   /*
155935aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
156035aab85fSBarry Smith     divisible by the blocksize
156135aab85fSBarry Smith   */
1562b6490206SBarry Smith   mbs        = M/bs;
156335aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
156435aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
156535aab85fSBarry Smith   else                  mbs++;
156635aab85fSBarry Smith   if (extra_rows) {
1567537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
156835aab85fSBarry Smith   }
1569b6490206SBarry Smith 
15702593348eSBarry Smith   /* read in row lengths */
157135aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
15720752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
157335aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
15742593348eSBarry Smith 
1575b6490206SBarry Smith   /* read in column indices */
157635aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
15770752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr);
157835aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1579b6490206SBarry Smith 
1580b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1581b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1582b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
158335aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
158435aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
158535aab85fSBarry Smith   masked = mask + mbs;
1586b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1587b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
158835aab85fSBarry Smith     nmask = 0;
1589b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1590b6490206SBarry Smith       kmax = rowlengths[rowcount];
1591b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
159235aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
159335aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1594b6490206SBarry Smith       }
1595b6490206SBarry Smith       rowcount++;
1596b6490206SBarry Smith     }
159735aab85fSBarry Smith     browlengths[i] += nmask;
159835aab85fSBarry Smith     /* zero out the mask elements we set */
159935aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1600b6490206SBarry Smith   }
1601b6490206SBarry Smith 
16022593348eSBarry Smith   /* create our matrix */
16033eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
16042593348eSBarry Smith   B = *A;
1605b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
16062593348eSBarry Smith 
16072593348eSBarry Smith   /* set matrix "i" values */
1608de6a44a3SBarry Smith   a->i[0] = 0;
1609b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1610b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1611b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
16122593348eSBarry Smith   }
1613b6490206SBarry Smith   a->nz         = 0;
1614b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
16152593348eSBarry Smith 
1616b6490206SBarry Smith   /* read in nonzero values */
161735aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
16180752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr);
161935aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1620b6490206SBarry Smith 
1621b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1622b6490206SBarry Smith   nzcount = 0; jcount = 0;
1623b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1624b6490206SBarry Smith     nzcountb = nzcount;
162535aab85fSBarry Smith     nmask    = 0;
1626b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1627b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1628b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
162935aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
163035aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1631b6490206SBarry Smith       }
1632b6490206SBarry Smith       rowcount++;
1633b6490206SBarry Smith     }
1634de6a44a3SBarry Smith     /* sort the masked values */
163577c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1636de6a44a3SBarry Smith 
1637b6490206SBarry Smith     /* set "j" values into matrix */
1638b6490206SBarry Smith     maskcount = 1;
163935aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
164035aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1641de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1642b6490206SBarry Smith     }
1643b6490206SBarry Smith     /* set "a" values into matrix */
1644de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1645b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1646b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1647b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1648de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
1649de6a44a3SBarry Smith         block  = mask[tmp] - 1;
1650de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
1651de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
1652b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1653b6490206SBarry Smith       }
1654b6490206SBarry Smith     }
165535aab85fSBarry Smith     /* zero out the mask elements we set */
165635aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1657b6490206SBarry Smith   }
1658a8c6a408SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1659b6490206SBarry Smith 
1660b6490206SBarry Smith   PetscFree(rowlengths);
1661b6490206SBarry Smith   PetscFree(browlengths);
1662b6490206SBarry Smith   PetscFree(aa);
1663b6490206SBarry Smith   PetscFree(jj);
1664b6490206SBarry Smith   PetscFree(mask);
1665b6490206SBarry Smith 
1666b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1667de6a44a3SBarry Smith 
16689c01be13SBarry Smith   ierr = MatView_Private(B); CHKERRQ(ierr);
16693a40ed3dSBarry Smith   PetscFunctionReturn(0);
16702593348eSBarry Smith }
16712593348eSBarry Smith 
16722593348eSBarry Smith 
16732593348eSBarry Smith 
1674