xref: /petsc/src/mat/impls/baij/seq/baij.c (revision bd648c37e6d565882dd6b09efc80b89e39bd7ace)
1*bd648c37SKris Buschelman /*$Id: baij.c,v 1.229 2001/06/21 23:47:10 buschelm Exp buschelm $*/
22593348eSBarry Smith 
32593348eSBarry Smith /*
4b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
52593348eSBarry Smith   matrix storage format.
62593348eSBarry Smith */
7a2ccead7SSatish Balay #include "petscsys.h"                     /*I "petscmat.h" I*/
870f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
91a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
101a6a6d98SBarry Smith #include "src/inline/spops.h"
113b547af2SSatish Balay 
1295d5f7c2SBarry Smith /*  UGLY, ugly, ugly
1395d5f7c2SBarry Smith    When MatScalar == Scalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does
1495d5f7c2SBarry Smith    not exist. Otherwise ..._MatScalar() takes matrix elements in single precision and
1595d5f7c2SBarry Smith    inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ()
1695d5f7c2SBarry Smith    converts the entries into single precision and then calls ..._MatScalar() to put them
1795d5f7c2SBarry Smith    into the single precision data structures.
1895d5f7c2SBarry Smith */
1995d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
20ca44d042SBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode);
2195d5f7c2SBarry Smith #else
2295d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ
2395d5f7c2SBarry Smith #endif
2495d5f7c2SBarry Smith 
252d61bbb3SSatish Balay #define CHUNKSIZE  10
26de6a44a3SBarry Smith 
27be5855fcSBarry Smith /*
28be5855fcSBarry Smith      Checks for missing diagonals
29be5855fcSBarry Smith */
304a2ae208SSatish Balay #undef __FUNCT__
314a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqBAIJ"
32c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A)
33be5855fcSBarry Smith {
34be5855fcSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
35883fce79SBarry Smith   int         *diag,*jj = a->j,i,ierr;
36be5855fcSBarry Smith 
37be5855fcSBarry Smith   PetscFunctionBegin;
38c4992f7dSBarry Smith   ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
39883fce79SBarry Smith   diag = a->diag;
400e8e8aceSBarry Smith   for (i=0; i<a->mbs; i++) {
41be5855fcSBarry Smith     if (jj[diag[i]] != i) {
4229bbc08cSBarry Smith       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
43be5855fcSBarry Smith     }
44be5855fcSBarry Smith   }
45be5855fcSBarry Smith   PetscFunctionReturn(0);
46be5855fcSBarry Smith }
47be5855fcSBarry Smith 
484a2ae208SSatish Balay #undef __FUNCT__
494a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqBAIJ"
50c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A)
51de6a44a3SBarry Smith {
52de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
5382502324SSatish Balay   int         i,j,*diag,m = a->mbs,ierr;
54de6a44a3SBarry Smith 
553a40ed3dSBarry Smith   PetscFunctionBegin;
56883fce79SBarry Smith   if (a->diag) PetscFunctionReturn(0);
57883fce79SBarry Smith 
58b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
59b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
607fc0212eSBarry Smith   for (i=0; i<m; i++) {
61ecc615b2SBarry Smith     diag[i] = a->i[i+1];
62de6a44a3SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
63de6a44a3SBarry Smith       if (a->j[j] == i) {
64de6a44a3SBarry Smith         diag[i] = j;
65de6a44a3SBarry Smith         break;
66de6a44a3SBarry Smith       }
67de6a44a3SBarry Smith     }
68de6a44a3SBarry Smith   }
69de6a44a3SBarry Smith   a->diag = diag;
703a40ed3dSBarry Smith   PetscFunctionReturn(0);
71de6a44a3SBarry Smith }
722593348eSBarry Smith 
732593348eSBarry Smith 
74ca44d042SBarry Smith EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
753b2fbd54SBarry Smith 
764a2ae208SSatish Balay #undef __FUNCT__
774a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqBAIJ"
780e6d2581SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
793b2fbd54SBarry Smith {
803b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
813b2fbd54SBarry Smith   int         ierr,n = a->mbs,i;
823b2fbd54SBarry Smith 
833a40ed3dSBarry Smith   PetscFunctionBegin;
843b2fbd54SBarry Smith   *nn = n;
853a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
863b2fbd54SBarry Smith   if (symmetric) {
873b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
883b2fbd54SBarry Smith   } else if (oshift == 1) {
893b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
90435da068SBarry Smith     int nz = a->i[n];
913b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
923b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]++;
933b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
943b2fbd54SBarry Smith   } else {
953b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
963b2fbd54SBarry Smith   }
973b2fbd54SBarry Smith 
983a40ed3dSBarry Smith   PetscFunctionReturn(0);
993b2fbd54SBarry Smith }
1003b2fbd54SBarry Smith 
1014a2ae208SSatish Balay #undef __FUNCT__
1024a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqBAIJ"
103435da068SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
1043b2fbd54SBarry Smith {
1053b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
106606d414cSSatish Balay   int         i,n = a->mbs,ierr;
1073b2fbd54SBarry Smith 
1083a40ed3dSBarry Smith   PetscFunctionBegin;
1093a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1103b2fbd54SBarry Smith   if (symmetric) {
111606d414cSSatish Balay     ierr = PetscFree(*ia);CHKERRQ(ierr);
112606d414cSSatish Balay     ierr = PetscFree(*ja);CHKERRQ(ierr);
113af5da2bfSBarry Smith   } else if (oshift == 1) {
114435da068SBarry Smith     int nz = a->i[n]-1;
1153b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
1163b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]--;
1173b2fbd54SBarry Smith   }
1183a40ed3dSBarry Smith   PetscFunctionReturn(0);
1193b2fbd54SBarry Smith }
1203b2fbd54SBarry Smith 
1214a2ae208SSatish Balay #undef __FUNCT__
1224a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqBAIJ"
1232d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs)
1242d61bbb3SSatish Balay {
1252d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data;
1262d61bbb3SSatish Balay 
1272d61bbb3SSatish Balay   PetscFunctionBegin;
1282d61bbb3SSatish Balay   *bs = baij->bs;
1292d61bbb3SSatish Balay   PetscFunctionReturn(0);
1302d61bbb3SSatish Balay }
1312d61bbb3SSatish Balay 
1322d61bbb3SSatish Balay 
1334a2ae208SSatish Balay #undef __FUNCT__
1344a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqBAIJ"
135e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
1362d61bbb3SSatish Balay {
1372d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
138e51c0b9cSSatish Balay   int         ierr;
1392d61bbb3SSatish Balay 
140433994e6SBarry Smith   PetscFunctionBegin;
141aa482453SBarry Smith #if defined(PETSC_USE_LOG)
142b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
1432d61bbb3SSatish Balay #endif
144606d414cSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
145606d414cSSatish Balay   if (!a->singlemalloc) {
146606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
147606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
148606d414cSSatish Balay   }
149c38d4ed2SBarry Smith   if (a->row) {
150c38d4ed2SBarry Smith     ierr = ISDestroy(a->row);CHKERRQ(ierr);
151c38d4ed2SBarry Smith   }
152c38d4ed2SBarry Smith   if (a->col) {
153c38d4ed2SBarry Smith     ierr = ISDestroy(a->col);CHKERRQ(ierr);
154c38d4ed2SBarry Smith   }
155606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
156606d414cSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
157606d414cSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
158606d414cSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
159606d414cSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
160e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
161606d414cSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
162563b5814SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
163563b5814SBarry Smith   if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);}
164563b5814SBarry Smith #endif
165606d414cSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
1662d61bbb3SSatish Balay   PetscFunctionReturn(0);
1672d61bbb3SSatish Balay }
1682d61bbb3SSatish Balay 
1694a2ae208SSatish Balay #undef __FUNCT__
1704a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqBAIJ"
1712d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1722d61bbb3SSatish Balay {
1732d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1742d61bbb3SSatish Balay 
1752d61bbb3SSatish Balay   PetscFunctionBegin;
176aa275fccSKris Buschelman   switch (op) {
177aa275fccSKris Buschelman   case MAT_ROW_ORIENTED:
178aa275fccSKris Buschelman     a->roworiented    = PETSC_TRUE;
179aa275fccSKris Buschelman     break;
180aa275fccSKris Buschelman   case MAT_COLUMN_ORIENTED:
181aa275fccSKris Buschelman     a->roworiented    = PETSC_FALSE;
182aa275fccSKris Buschelman     break;
183aa275fccSKris Buschelman   case MAT_COLUMNS_SORTED:
184aa275fccSKris Buschelman     a->sorted         = PETSC_TRUE;
185aa275fccSKris Buschelman     break;
186aa275fccSKris Buschelman   case MAT_COLUMNS_UNSORTED:
187aa275fccSKris Buschelman     a->sorted         = PETSC_FALSE;
188aa275fccSKris Buschelman     break;
189aa275fccSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
190aa275fccSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
191aa275fccSKris Buschelman     break;
192aa275fccSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
193aa275fccSKris Buschelman     a->nonew          = 1;
194aa275fccSKris Buschelman     break;
195aa275fccSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
196aa275fccSKris Buschelman     a->nonew          = -1;
197aa275fccSKris Buschelman     break;
198aa275fccSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
199aa275fccSKris Buschelman     a->nonew          = -2;
200aa275fccSKris Buschelman     break;
201aa275fccSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
202aa275fccSKris Buschelman     a->nonew          = 0;
203aa275fccSKris Buschelman     break;
204aa275fccSKris Buschelman   case MAT_ROWS_SORTED:
205aa275fccSKris Buschelman   case MAT_ROWS_UNSORTED:
206aa275fccSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
207aa275fccSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
208aa275fccSKris Buschelman   case MAT_USE_HASH_TABLE:
209b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
210aa275fccSKris Buschelman     break;
211aa275fccSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
21229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
213*bd648c37SKris Buschelman     break;
214*bd648c37SKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
215*bd648c37SKris Buschelman #if defined(PETSC_USE_MAT_SINGLE)
216*bd648c37SKris Buschelman     if (a->bs==4) {
217*bd648c37SKris Buschelman       PetscTruth sse_enabled;
218*bd648c37SKris Buschelman       int ierr;
219*bd648c37SKris Buschelman       ierr = PetscSSEIsEnabled(&sse_enabled);CHKERRQ(ierr);
220*bd648c37SKris Buschelman       if (sse_enabled) {
221*bd648c37SKris Buschelman         singleprecisionsolves = PETSC_TRUE;
222*bd648c37SKris Buschelman       } else {
223*bd648c37SKris Buschelman         PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n");
224*bd648c37SKris Buschelman       }
225*bd648c37SKris Buschelman     } else {
226*bd648c37SKris Buschelman       PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n");
227*bd648c37SKris Buschelman     }
228*bd648c37SKris Buschelman #else
229*bd648c37SKris Buschelman     PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n");
230*bd648c37SKris Buschelman #endif
231*bd648c37SKris Buschelman     break;
232aa275fccSKris Buschelman   default:
23329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
2342d61bbb3SSatish Balay   }
2352d61bbb3SSatish Balay   PetscFunctionReturn(0);
2362d61bbb3SSatish Balay }
2372d61bbb3SSatish Balay 
2384a2ae208SSatish Balay #undef __FUNCT__
2394a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_SeqBAIJ"
2402d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
2412d61bbb3SSatish Balay {
2422d61bbb3SSatish Balay   PetscFunctionBegin;
2434c49b128SBarry Smith   if (m) *m = 0;
244273d9f13SBarry Smith   if (n) *n = A->m;
2452d61bbb3SSatish Balay   PetscFunctionReturn(0);
2462d61bbb3SSatish Balay }
2472d61bbb3SSatish Balay 
2484a2ae208SSatish Balay #undef __FUNCT__
2494a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ"
2502d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2512d61bbb3SSatish Balay {
2522d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
25382502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr;
2543f1db9ecSBarry Smith   MatScalar    *aa,*aa_i;
2553f1db9ecSBarry Smith   Scalar       *v_i;
2562d61bbb3SSatish Balay 
2572d61bbb3SSatish Balay   PetscFunctionBegin;
2582d61bbb3SSatish Balay   bs  = a->bs;
2592d61bbb3SSatish Balay   ai  = a->i;
2602d61bbb3SSatish Balay   aj  = a->j;
2612d61bbb3SSatish Balay   aa  = a->a;
2622d61bbb3SSatish Balay   bs2 = a->bs2;
2632d61bbb3SSatish Balay 
264273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
2652d61bbb3SSatish Balay 
2662d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2672d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2682d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2692d61bbb3SSatish Balay   *nz = bs*M;
2702d61bbb3SSatish Balay 
2712d61bbb3SSatish Balay   if (v) {
2722d61bbb3SSatish Balay     *v = 0;
2732d61bbb3SSatish Balay     if (*nz) {
274b0a32e0cSBarry Smith       ierr = PetscMalloc((*nz)*sizeof(Scalar),v);CHKERRQ(ierr);
2752d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
2762d61bbb3SSatish Balay         v_i  = *v + i*bs;
2772d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2782d61bbb3SSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
2792d61bbb3SSatish Balay       }
2802d61bbb3SSatish Balay     }
2812d61bbb3SSatish Balay   }
2822d61bbb3SSatish Balay 
2832d61bbb3SSatish Balay   if (idx) {
2842d61bbb3SSatish Balay     *idx = 0;
2852d61bbb3SSatish Balay     if (*nz) {
286b0a32e0cSBarry Smith       ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr);
2872d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
2882d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2892d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2902d61bbb3SSatish Balay         for (j=0; j<bs; j++) {idx_i[j] = itmp++;}
2912d61bbb3SSatish Balay       }
2922d61bbb3SSatish Balay     }
2932d61bbb3SSatish Balay   }
2942d61bbb3SSatish Balay   PetscFunctionReturn(0);
2952d61bbb3SSatish Balay }
2962d61bbb3SSatish Balay 
2974a2ae208SSatish Balay #undef __FUNCT__
2984a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ"
2992d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
3002d61bbb3SSatish Balay {
301606d414cSSatish Balay   int ierr;
302606d414cSSatish Balay 
3032d61bbb3SSatish Balay   PetscFunctionBegin;
304606d414cSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
305606d414cSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
3062d61bbb3SSatish Balay   PetscFunctionReturn(0);
3072d61bbb3SSatish Balay }
3082d61bbb3SSatish Balay 
3094a2ae208SSatish Balay #undef __FUNCT__
3104a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ"
3112d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
3122d61bbb3SSatish Balay {
3132d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
3142d61bbb3SSatish Balay   Mat         C;
3152d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
316273d9f13SBarry Smith   int         *rows,*cols,bs2=a->bs2;
317375fe846SBarry Smith   Scalar      *array;
3182d61bbb3SSatish Balay 
3192d61bbb3SSatish Balay   PetscFunctionBegin;
32029bbc08cSBarry Smith   if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place");
321b0a32e0cSBarry Smith   ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr);
322549d3d68SSatish Balay   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
3232d61bbb3SSatish Balay 
324375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
325b0a32e0cSBarry Smith   ierr = PetscMalloc(a->bs2*a->nz*sizeof(Scalar),&array);CHKERRQ(ierr);
326375fe846SBarry Smith   for (i=0; i<a->bs2*a->nz; i++) array[i] = (Scalar)a->a[i];
327375fe846SBarry Smith #else
3283eda8832SBarry Smith   array = a->a;
329375fe846SBarry Smith #endif
330375fe846SBarry Smith 
3312d61bbb3SSatish Balay   for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1;
332273d9f13SBarry Smith   ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
333606d414cSSatish Balay   ierr = PetscFree(col);CHKERRQ(ierr);
334b0a32e0cSBarry Smith   ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr);
3352d61bbb3SSatish Balay   cols = rows + bs;
3362d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3372d61bbb3SSatish Balay     cols[0] = i*bs;
3382d61bbb3SSatish Balay     for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1;
3392d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
3402d61bbb3SSatish Balay     for (j=0; j<len; j++) {
3412d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
3422d61bbb3SSatish Balay       for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1;
3432d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
3442d61bbb3SSatish Balay       array += bs2;
3452d61bbb3SSatish Balay     }
3462d61bbb3SSatish Balay   }
347606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
348375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
349375fe846SBarry Smith   ierr = PetscFree(array);
350375fe846SBarry Smith #endif
3512d61bbb3SSatish Balay 
3522d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3532d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3542d61bbb3SSatish Balay 
355c4992f7dSBarry Smith   if (B) {
3562d61bbb3SSatish Balay     *B = C;
3572d61bbb3SSatish Balay   } else {
358273d9f13SBarry Smith     ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
3592d61bbb3SSatish Balay   }
3602d61bbb3SSatish Balay   PetscFunctionReturn(0);
3612d61bbb3SSatish Balay }
3622d61bbb3SSatish Balay 
3634a2ae208SSatish Balay #undef __FUNCT__
3644a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary"
365b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer)
3662593348eSBarry Smith {
367b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
3689df24120SSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
369b6490206SBarry Smith   Scalar      *aa;
370ce6f0cecSBarry Smith   FILE        *file;
3712593348eSBarry Smith 
3723a40ed3dSBarry Smith   PetscFunctionBegin;
373b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
374b0a32e0cSBarry Smith   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
3752593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3763b2fbd54SBarry Smith 
377273d9f13SBarry Smith   col_lens[1] = A->m;
378273d9f13SBarry Smith   col_lens[2] = A->n;
3797e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3802593348eSBarry Smith 
3812593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
382b6490206SBarry Smith   count = 0;
383b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
384b6490206SBarry Smith     for (j=0; j<bs; j++) {
385b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
386b6490206SBarry Smith     }
3872593348eSBarry Smith   }
388273d9f13SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
389606d414cSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
3902593348eSBarry Smith 
3912593348eSBarry Smith   /* store column indices (zero start index) */
392b0a32e0cSBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
393b6490206SBarry Smith   count = 0;
394b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
395b6490206SBarry Smith     for (j=0; j<bs; j++) {
396b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
397b6490206SBarry Smith         for (l=0; l<bs; l++) {
398b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3992593348eSBarry Smith         }
4002593348eSBarry Smith       }
401b6490206SBarry Smith     }
402b6490206SBarry Smith   }
4030752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
404606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
4052593348eSBarry Smith 
4062593348eSBarry Smith   /* store nonzero values */
407b0a32e0cSBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(Scalar),&aa);CHKERRQ(ierr);
408b6490206SBarry Smith   count = 0;
409b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
410b6490206SBarry Smith     for (j=0; j<bs; j++) {
411b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
412b6490206SBarry Smith         for (l=0; l<bs; l++) {
4137e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
414b6490206SBarry Smith         }
415b6490206SBarry Smith       }
416b6490206SBarry Smith     }
417b6490206SBarry Smith   }
4180752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
419606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
420ce6f0cecSBarry Smith 
421b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
422ce6f0cecSBarry Smith   if (file) {
423ce6f0cecSBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
424ce6f0cecSBarry Smith   }
4253a40ed3dSBarry Smith   PetscFunctionReturn(0);
4262593348eSBarry Smith }
4272593348eSBarry Smith 
4284a2ae208SSatish Balay #undef __FUNCT__
4294a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII"
430b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer)
4312593348eSBarry Smith {
432b6490206SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
433fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
434f3ef73ceSBarry Smith   PetscViewerFormat format;
4352593348eSBarry Smith 
4363a40ed3dSBarry Smith   PetscFunctionBegin;
437b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
438fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
439b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
440fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
44129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
442fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
443b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
44444cd7ae7SLois Curfman McInnes     for (i=0; i<a->mbs; i++) {
44544cd7ae7SLois Curfman McInnes       for (j=0; j<bs; j++) {
446b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
44744cd7ae7SLois Curfman McInnes         for (k=a->i[i]; k<a->i[i+1]; k++) {
44844cd7ae7SLois Curfman McInnes           for (l=0; l<bs; l++) {
449aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4500e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
451b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %gi",bs*a->j[k]+l,
4520e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4530e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
454b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %gi",bs*a->j[k]+l,
4550e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4560e6d2581SBarry Smith             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
457b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4580ef38995SBarry Smith             }
45944cd7ae7SLois Curfman McInnes #else
4600ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
461b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
4620ef38995SBarry Smith             }
46344cd7ae7SLois Curfman McInnes #endif
46444cd7ae7SLois Curfman McInnes           }
46544cd7ae7SLois Curfman McInnes         }
466b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
46744cd7ae7SLois Curfman McInnes       }
46844cd7ae7SLois Curfman McInnes     }
469b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
4700ef38995SBarry Smith   } else {
471b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
472b6490206SBarry Smith     for (i=0; i<a->mbs; i++) {
473b6490206SBarry Smith       for (j=0; j<bs; j++) {
474b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
475b6490206SBarry Smith         for (k=a->i[i]; k<a->i[i+1]; k++) {
476b6490206SBarry Smith           for (l=0; l<bs; l++) {
477aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4780e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
479b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
4800e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4810e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
482b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
4830e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4840ef38995SBarry Smith             } else {
485b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
48688685aaeSLois Curfman McInnes             }
48788685aaeSLois Curfman McInnes #else
488b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
48988685aaeSLois Curfman McInnes #endif
4902593348eSBarry Smith           }
4912593348eSBarry Smith         }
492b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
4932593348eSBarry Smith       }
4942593348eSBarry Smith     }
495b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
496b6490206SBarry Smith   }
497b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
4983a40ed3dSBarry Smith   PetscFunctionReturn(0);
4992593348eSBarry Smith }
5002593348eSBarry Smith 
5014a2ae208SSatish Balay #undef __FUNCT__
5024a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom"
503b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
5043270192aSSatish Balay {
50577ed5343SBarry Smith   Mat          A = (Mat) Aa;
5063270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ*)A->data;
507b65db4caSBarry Smith   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2;
5080e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
5093f1db9ecSBarry Smith   MatScalar    *aa;
510b0a32e0cSBarry Smith   PetscViewer  viewer;
5113270192aSSatish Balay 
5123a40ed3dSBarry Smith   PetscFunctionBegin;
5133270192aSSatish Balay 
514b65db4caSBarry Smith   /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/
51577ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
51677ed5343SBarry Smith 
517b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
51877ed5343SBarry Smith 
5193270192aSSatish Balay   /* loop over matrix elements drawing boxes */
520b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
5213270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5223270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
523273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5243270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5253270192aSSatish Balay       aa = a->a + j*bs2;
5263270192aSSatish Balay       for (k=0; k<bs; k++) {
5273270192aSSatish Balay         for (l=0; l<bs; l++) {
5280e6d2581SBarry Smith           if (PetscRealPart(*aa++) >=  0.) continue;
529b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5303270192aSSatish Balay         }
5313270192aSSatish Balay       }
5323270192aSSatish Balay     }
5333270192aSSatish Balay   }
534b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
5353270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5363270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
537273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5383270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5393270192aSSatish Balay       aa = a->a + j*bs2;
5403270192aSSatish Balay       for (k=0; k<bs; k++) {
5413270192aSSatish Balay         for (l=0; l<bs; l++) {
5420e6d2581SBarry Smith           if (PetscRealPart(*aa++) != 0.) continue;
543b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5443270192aSSatish Balay         }
5453270192aSSatish Balay       }
5463270192aSSatish Balay     }
5473270192aSSatish Balay   }
5483270192aSSatish Balay 
549b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
5503270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5513270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
552273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5533270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5543270192aSSatish Balay       aa = a->a + j*bs2;
5553270192aSSatish Balay       for (k=0; k<bs; k++) {
5563270192aSSatish Balay         for (l=0; l<bs; l++) {
5570e6d2581SBarry Smith           if (PetscRealPart(*aa++) <= 0.) continue;
558b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5593270192aSSatish Balay         }
5603270192aSSatish Balay       }
5613270192aSSatish Balay     }
5623270192aSSatish Balay   }
56377ed5343SBarry Smith   PetscFunctionReturn(0);
56477ed5343SBarry Smith }
5653270192aSSatish Balay 
5664a2ae208SSatish Balay #undef __FUNCT__
5674a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw"
568b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer)
56977ed5343SBarry Smith {
5707dce120fSSatish Balay   int          ierr;
5710e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,w,h;
572b0a32e0cSBarry Smith   PetscDraw    draw;
57377ed5343SBarry Smith   PetscTruth   isnull;
5743270192aSSatish Balay 
57577ed5343SBarry Smith   PetscFunctionBegin;
57677ed5343SBarry Smith 
577b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
578b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
57977ed5343SBarry Smith 
58077ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
581273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
58277ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
583b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
584b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
58577ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5863a40ed3dSBarry Smith   PetscFunctionReturn(0);
5873270192aSSatish Balay }
5883270192aSSatish Balay 
5894a2ae208SSatish Balay #undef __FUNCT__
5904a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ"
591b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer)
5922593348eSBarry Smith {
59319bcc07fSBarry Smith   int        ierr;
5946831982aSBarry Smith   PetscTruth issocket,isascii,isbinary,isdraw;
5952593348eSBarry Smith 
5963a40ed3dSBarry Smith   PetscFunctionBegin;
597b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
598b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
599fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
600fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6010f5bd95cSBarry Smith   if (issocket) {
60229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
6030f5bd95cSBarry Smith   } else if (isascii){
6043a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
6050f5bd95cSBarry Smith   } else if (isbinary) {
6063a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
6070f5bd95cSBarry Smith   } else if (isdraw) {
6083a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
6095cd90555SBarry Smith   } else {
61029bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
6112593348eSBarry Smith   }
6123a40ed3dSBarry Smith   PetscFunctionReturn(0);
6132593348eSBarry Smith }
614b6490206SBarry Smith 
615cd0e1443SSatish Balay 
6164a2ae208SSatish Balay #undef __FUNCT__
6174a2ae208SSatish Balay #define __FUNCT__ "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;
6243f1db9ecSBarry Smith   MatScalar  *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;
62929bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
630273d9f13SBarry Smith     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"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 */
63429bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
635273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"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 
66195d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
6624a2ae208SSatish Balay #undef __FUNCT__
6634a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
66495d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv)
66595d5f7c2SBarry Smith {
666563b5814SBarry Smith   Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data;
667563b5814SBarry Smith   int         ierr,i,N = m*n*b->bs2;
668563b5814SBarry Smith   MatScalar   *vsingle;
66995d5f7c2SBarry Smith 
67095d5f7c2SBarry Smith   PetscFunctionBegin;
671563b5814SBarry Smith   if (N > b->setvalueslen) {
672563b5814SBarry Smith     if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);}
673b0a32e0cSBarry Smith     ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr);
674563b5814SBarry Smith     b->setvalueslen  = N;
675563b5814SBarry Smith   }
676563b5814SBarry Smith   vsingle = b->setvaluescopy;
67795d5f7c2SBarry Smith   for (i=0; i<N; i++) {
67895d5f7c2SBarry Smith     vsingle[i] = v[i];
67995d5f7c2SBarry Smith   }
68095d5f7c2SBarry Smith   ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr);
68195d5f7c2SBarry Smith   PetscFunctionReturn(0);
68295d5f7c2SBarry Smith }
68395d5f7c2SBarry Smith #endif
68495d5f7c2SBarry Smith 
6852d61bbb3SSatish Balay 
6864a2ae208SSatish Balay #undef __FUNCT__
6874a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
68895d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is)
68992c4ed94SBarry Smith {
69092c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6918a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
692273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
693549d3d68SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
694273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
69595d5f7c2SBarry Smith   MatScalar   *value = v,*ap,*aa = a->a,*bap;
69692c4ed94SBarry Smith 
6973a40ed3dSBarry Smith   PetscFunctionBegin;
6980e324ae4SSatish Balay   if (roworiented) {
6990e324ae4SSatish Balay     stepval = (n-1)*bs;
7000e324ae4SSatish Balay   } else {
7010e324ae4SSatish Balay     stepval = (m-1)*bs;
7020e324ae4SSatish Balay   }
70392c4ed94SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
70492c4ed94SBarry Smith     row  = im[k];
7055ef9f2a5SBarry Smith     if (row < 0) continue;
706aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
70729bbc08cSBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
70892c4ed94SBarry Smith #endif
70992c4ed94SBarry Smith     rp   = aj + ai[row];
71092c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
71192c4ed94SBarry Smith     rmax = imax[row];
71292c4ed94SBarry Smith     nrow = ailen[row];
71392c4ed94SBarry Smith     low  = 0;
71492c4ed94SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
7155ef9f2a5SBarry Smith       if (in[l] < 0) continue;
716aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
71729bbc08cSBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
71892c4ed94SBarry Smith #endif
71992c4ed94SBarry Smith       col = in[l];
72092c4ed94SBarry Smith       if (roworiented) {
7210e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
7220e324ae4SSatish Balay       } else {
7230e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
72492c4ed94SBarry Smith       }
72592c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
72692c4ed94SBarry Smith       while (high-low > 7) {
72792c4ed94SBarry Smith         t = (low+high)/2;
72892c4ed94SBarry Smith         if (rp[t] > col) high = t;
72992c4ed94SBarry Smith         else             low  = t;
73092c4ed94SBarry Smith       }
73192c4ed94SBarry Smith       for (i=low; i<high; i++) {
73292c4ed94SBarry Smith         if (rp[i] > col) break;
73392c4ed94SBarry Smith         if (rp[i] == col) {
7348a84c255SSatish Balay           bap  = ap +  bs2*i;
7350e324ae4SSatish Balay           if (roworiented) {
7368a84c255SSatish Balay             if (is == ADD_VALUES) {
737dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
738dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
7398a84c255SSatish Balay                   bap[jj] += *value++;
740dd9472c6SBarry Smith                 }
741dd9472c6SBarry Smith               }
7420e324ae4SSatish Balay             } else {
743dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
744dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
7450e324ae4SSatish Balay                   bap[jj] = *value++;
7468a84c255SSatish Balay                 }
747dd9472c6SBarry Smith               }
748dd9472c6SBarry Smith             }
7490e324ae4SSatish Balay           } else {
7500e324ae4SSatish Balay             if (is == ADD_VALUES) {
751dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
752dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
7530e324ae4SSatish Balay                   *bap++ += *value++;
754dd9472c6SBarry Smith                 }
755dd9472c6SBarry Smith               }
7560e324ae4SSatish Balay             } else {
757dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
758dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
7590e324ae4SSatish Balay                   *bap++  = *value++;
7600e324ae4SSatish Balay                 }
7618a84c255SSatish Balay               }
762dd9472c6SBarry Smith             }
763dd9472c6SBarry Smith           }
764f1241b54SBarry Smith           goto noinsert2;
76592c4ed94SBarry Smith         }
76692c4ed94SBarry Smith       }
76789280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
76829bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
76992c4ed94SBarry Smith       if (nrow >= rmax) {
77092c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
77192c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
7723f1db9ecSBarry Smith         MatScalar *new_a;
77392c4ed94SBarry Smith 
77429bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
77589280ab3SLois Curfman McInnes 
77692c4ed94SBarry Smith         /* malloc new storage space */
7773f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
778b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
77992c4ed94SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
78092c4ed94SBarry Smith         new_i   = new_j + new_nz;
78192c4ed94SBarry Smith 
78292c4ed94SBarry Smith         /* copy over old data into new slots */
78392c4ed94SBarry Smith         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
78492c4ed94SBarry Smith         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
785549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
78692c4ed94SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
787549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
788549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
789549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
790549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
79192c4ed94SBarry Smith         /* free up old matrix storage */
792606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
793606d414cSSatish Balay         if (!a->singlemalloc) {
794606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
795606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
796606d414cSSatish Balay         }
79792c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
798c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
79992c4ed94SBarry Smith 
80092c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
80192c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
802b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
80392c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
80492c4ed94SBarry Smith         a->reallocs++;
80592c4ed94SBarry Smith         a->nz++;
80692c4ed94SBarry Smith       }
80792c4ed94SBarry Smith       N = nrow++ - 1;
80892c4ed94SBarry Smith       /* shift up all the later entries in this row */
80992c4ed94SBarry Smith       for (ii=N; ii>=i; ii--) {
81092c4ed94SBarry Smith         rp[ii+1] = rp[ii];
811549d3d68SSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
81292c4ed94SBarry Smith       }
813549d3d68SSatish Balay       if (N >= i) {
814549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
815549d3d68SSatish Balay       }
81692c4ed94SBarry Smith       rp[i] = col;
8178a84c255SSatish Balay       bap   = ap +  bs2*i;
8180e324ae4SSatish Balay       if (roworiented) {
819dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
820dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs) {
8210e324ae4SSatish Balay             bap[jj] = *value++;
822dd9472c6SBarry Smith           }
823dd9472c6SBarry Smith         }
8240e324ae4SSatish Balay       } else {
825dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
826dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++) {
8270e324ae4SSatish Balay             *bap++  = *value++;
8280e324ae4SSatish Balay           }
829dd9472c6SBarry Smith         }
830dd9472c6SBarry Smith       }
831f1241b54SBarry Smith       noinsert2:;
83292c4ed94SBarry Smith       low = i;
83392c4ed94SBarry Smith     }
83492c4ed94SBarry Smith     ailen[row] = nrow;
83592c4ed94SBarry Smith   }
8363a40ed3dSBarry Smith   PetscFunctionReturn(0);
83792c4ed94SBarry Smith }
83892c4ed94SBarry Smith 
839f2501298SSatish Balay 
8404a2ae208SSatish Balay #undef __FUNCT__
8414a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ"
8428f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
843584200bdSSatish Balay {
844584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
845584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
846273d9f13SBarry Smith   int        m = A->m,*ip,N,*ailen = a->ilen;
847549d3d68SSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
8483f1db9ecSBarry Smith   MatScalar  *aa = a->a,*ap;
849584200bdSSatish Balay 
8503a40ed3dSBarry Smith   PetscFunctionBegin;
8513a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
852584200bdSSatish Balay 
85343ee02c3SBarry Smith   if (m) rmax = ailen[0];
854584200bdSSatish Balay   for (i=1; i<mbs; i++) {
855584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
856584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
857d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
858584200bdSSatish Balay     if (fshift) {
859a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
860584200bdSSatish Balay       N = ailen[i];
861584200bdSSatish Balay       for (j=0; j<N; j++) {
862584200bdSSatish Balay         ip[j-fshift] = ip[j];
863549d3d68SSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
864584200bdSSatish Balay       }
865584200bdSSatish Balay     }
866584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
867584200bdSSatish Balay   }
868584200bdSSatish Balay   if (mbs) {
869584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
870584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
871584200bdSSatish Balay   }
872584200bdSSatish Balay   /* reset ilen and imax for each row */
873584200bdSSatish Balay   for (i=0; i<mbs; i++) {
874584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
875584200bdSSatish Balay   }
876a7c10996SSatish Balay   a->nz = ai[mbs];
877584200bdSSatish Balay 
878584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
879584200bdSSatish Balay   if (fshift && a->diag) {
880606d414cSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
881b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-(m+1)*sizeof(int));
882584200bdSSatish Balay     a->diag = 0;
883584200bdSSatish Balay   }
884b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
885273d9f13SBarry Smith            m,A->n,a->bs,fshift*bs2,a->nz*bs2);
886b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
887584200bdSSatish Balay            a->reallocs);
888b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
889e2f3b5e9SSatish Balay   a->reallocs          = 0;
8900e6d2581SBarry Smith   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
8914e220ebcSLois Curfman McInnes 
8923a40ed3dSBarry Smith   PetscFunctionReturn(0);
893584200bdSSatish Balay }
894584200bdSSatish Balay 
8955a838052SSatish Balay 
896bea157c4SSatish Balay 
897bea157c4SSatish Balay /*
898bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
899bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
900bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
901bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
902bea157c4SSatish Balay */
9034a2ae208SSatish Balay #undef __FUNCT__
9044a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks"
905bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
906d9b7c43dSSatish Balay {
907bea157c4SSatish Balay   int        i,j,k,row;
908bea157c4SSatish Balay   PetscTruth flg;
9093a40ed3dSBarry Smith 
910433994e6SBarry Smith   PetscFunctionBegin;
911bea157c4SSatish Balay   for (i=0,j=0; i<n; j++) {
912bea157c4SSatish Balay     row = idx[i];
913bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
914bea157c4SSatish Balay       sizes[j] = 1;
915bea157c4SSatish Balay       i++;
916e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
917bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
918bea157c4SSatish Balay       i++;
919bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
920bea157c4SSatish Balay       flg = PETSC_TRUE;
921bea157c4SSatish Balay       for (k=1; k<bs; k++) {
922bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
923bea157c4SSatish Balay           flg = PETSC_FALSE;
924bea157c4SSatish Balay           break;
925d9b7c43dSSatish Balay         }
926bea157c4SSatish Balay       }
927bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
928bea157c4SSatish Balay         sizes[j] = bs;
929bea157c4SSatish Balay         i+= bs;
930bea157c4SSatish Balay       } else {
931bea157c4SSatish Balay         sizes[j] = 1;
932bea157c4SSatish Balay         i++;
933bea157c4SSatish Balay       }
934bea157c4SSatish Balay     }
935bea157c4SSatish Balay   }
936bea157c4SSatish Balay   *bs_max = j;
9373a40ed3dSBarry Smith   PetscFunctionReturn(0);
938d9b7c43dSSatish Balay }
939d9b7c43dSSatish Balay 
9404a2ae208SSatish Balay #undef __FUNCT__
9414a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ"
9428f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is,Scalar *diag)
943d9b7c43dSSatish Balay {
944d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
945b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
946bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
9473f1db9ecSBarry Smith   Scalar      zero = 0.0;
9483f1db9ecSBarry Smith   MatScalar   *aa;
949d9b7c43dSSatish Balay 
9503a40ed3dSBarry Smith   PetscFunctionBegin;
951d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
952b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr);
953d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
954d9b7c43dSSatish Balay 
955bea157c4SSatish Balay   /* allocate memory for rows,sizes */
956b0a32e0cSBarry Smith   ierr  = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
957bea157c4SSatish Balay   sizes = rows + is_n;
958bea157c4SSatish Balay 
959563b5814SBarry Smith   /* copy IS values to rows, and sort them */
960bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
961bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
962dffd3267SBarry Smith   if (baij->keepzeroedrows) {
963dffd3267SBarry Smith     for (i=0; i<is_n; i++) { sizes[i] = 1; }
964dffd3267SBarry Smith     bs_max = is_n;
965dffd3267SBarry Smith   } else {
966bea157c4SSatish Balay     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
967dffd3267SBarry Smith   }
968b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
969bea157c4SSatish Balay 
970bea157c4SSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
971bea157c4SSatish Balay     row   = rows[j];
972273d9f13SBarry Smith     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
973bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
974bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
975dffd3267SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
976bea157c4SSatish Balay       if (diag) {
977bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
978bea157c4SSatish Balay           baij->ilen[row/bs]       = 1;
979bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
980bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
981a07cd24cSSatish Balay         }
982563b5814SBarry Smith         /* Now insert all the diagonal values for this bs */
983bea157c4SSatish Balay         for (k=0; k<bs; k++) {
984bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
985bea157c4SSatish Balay         }
986bea157c4SSatish Balay       } else { /* (!diag) */
987bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
988bea157c4SSatish Balay       } /* end (!diag) */
989bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
990aa482453SBarry Smith #if defined (PETSC_USE_DEBUG)
99129bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
992bea157c4SSatish Balay #endif
993bea157c4SSatish Balay       for (k=0; k<count; k++) {
994d9b7c43dSSatish Balay         aa[0] =  zero;
995d9b7c43dSSatish Balay         aa    += bs;
996d9b7c43dSSatish Balay       }
997d9b7c43dSSatish Balay       if (diag) {
998f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
999d9b7c43dSSatish Balay       }
1000d9b7c43dSSatish Balay     }
1001bea157c4SSatish Balay   }
1002bea157c4SSatish Balay 
1003606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
10049a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10053a40ed3dSBarry Smith   PetscFunctionReturn(0);
1006d9b7c43dSSatish Balay }
10071c351548SSatish Balay 
10084a2ae208SSatish Balay #undef __FUNCT__
10094a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ"
10102d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
10112d61bbb3SSatish Balay {
10122d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
10132d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
1014273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
10152d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1016549d3d68SSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
1017273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
10183f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
10192d61bbb3SSatish Balay 
10202d61bbb3SSatish Balay   PetscFunctionBegin;
10212d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
10222d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
10235ef9f2a5SBarry Smith     if (row < 0) continue;
1024aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1025273d9f13SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
10262d61bbb3SSatish Balay #endif
10272d61bbb3SSatish Balay     rp   = aj + ai[brow];
10282d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
10292d61bbb3SSatish Balay     rmax = imax[brow];
10302d61bbb3SSatish Balay     nrow = ailen[brow];
10312d61bbb3SSatish Balay     low  = 0;
10322d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
10335ef9f2a5SBarry Smith       if (in[l] < 0) continue;
1034aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1035273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n);
10362d61bbb3SSatish Balay #endif
10372d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
10382d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
10392d61bbb3SSatish Balay       if (roworiented) {
10405ef9f2a5SBarry Smith         value = v[l + k*n];
10412d61bbb3SSatish Balay       } else {
10422d61bbb3SSatish Balay         value = v[k + l*m];
10432d61bbb3SSatish Balay       }
10442d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
10452d61bbb3SSatish Balay       while (high-low > 7) {
10462d61bbb3SSatish Balay         t = (low+high)/2;
10472d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
10482d61bbb3SSatish Balay         else              low  = t;
10492d61bbb3SSatish Balay       }
10502d61bbb3SSatish Balay       for (i=low; i<high; i++) {
10512d61bbb3SSatish Balay         if (rp[i] > bcol) break;
10522d61bbb3SSatish Balay         if (rp[i] == bcol) {
10532d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
10542d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
10552d61bbb3SSatish Balay           else                  *bap  = value;
10562d61bbb3SSatish Balay           goto noinsert1;
10572d61bbb3SSatish Balay         }
10582d61bbb3SSatish Balay       }
10592d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
106029bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
10612d61bbb3SSatish Balay       if (nrow >= rmax) {
10622d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
10632d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
10643f1db9ecSBarry Smith         MatScalar *new_a;
10652d61bbb3SSatish Balay 
106629bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
10672d61bbb3SSatish Balay 
10682d61bbb3SSatish Balay         /* Malloc new storage space */
10693f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
1070b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
10712d61bbb3SSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
10722d61bbb3SSatish Balay         new_i   = new_j + new_nz;
10732d61bbb3SSatish Balay 
10742d61bbb3SSatish Balay         /* copy over old data into new slots */
10752d61bbb3SSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
10762d61bbb3SSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1077549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
10782d61bbb3SSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1079549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1080549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1081549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1082549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
10832d61bbb3SSatish Balay         /* free up old matrix storage */
1084606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
1085606d414cSSatish Balay         if (!a->singlemalloc) {
1086606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
1087606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
1088606d414cSSatish Balay         }
10892d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1090c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
10912d61bbb3SSatish Balay 
10922d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10932d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1094b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
10952d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10962d61bbb3SSatish Balay         a->reallocs++;
10972d61bbb3SSatish Balay         a->nz++;
10982d61bbb3SSatish Balay       }
10992d61bbb3SSatish Balay       N = nrow++ - 1;
11002d61bbb3SSatish Balay       /* shift up all the later entries in this row */
11012d61bbb3SSatish Balay       for (ii=N; ii>=i; ii--) {
11022d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
1103549d3d68SSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
11042d61bbb3SSatish Balay       }
1105549d3d68SSatish Balay       if (N>=i) {
1106549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1107549d3d68SSatish Balay       }
11082d61bbb3SSatish Balay       rp[i]                      = bcol;
11092d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
11102d61bbb3SSatish Balay       noinsert1:;
11112d61bbb3SSatish Balay       low = i;
11122d61bbb3SSatish Balay     }
11132d61bbb3SSatish Balay     ailen[brow] = nrow;
11142d61bbb3SSatish Balay   }
11152d61bbb3SSatish Balay   PetscFunctionReturn(0);
11162d61bbb3SSatish Balay }
11172d61bbb3SSatish Balay 
1118b9b97703SBarry Smith EXTERN int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatLUInfo*,Mat*);
1119b9b97703SBarry Smith EXTERN int MatLUFactor_SeqBAIJ(Mat,IS,IS,MatLUInfo*);
1120ca44d042SBarry Smith EXTERN int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
1121ca44d042SBarry Smith EXTERN int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
1122ca44d042SBarry Smith EXTERN int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
1123ca44d042SBarry Smith EXTERN int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
1124ca44d042SBarry Smith EXTERN int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
1125ca44d042SBarry Smith EXTERN int MatScale_SeqBAIJ(Scalar*,Mat);
1126ca44d042SBarry Smith EXTERN int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
1127ca44d042SBarry Smith EXTERN int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
1128ca44d042SBarry Smith EXTERN int MatGetDiagonal_SeqBAIJ(Mat,Vec);
1129ca44d042SBarry Smith EXTERN int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
1130ca44d042SBarry Smith EXTERN int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
1131ca44d042SBarry Smith EXTERN int MatZeroEntries_SeqBAIJ(Mat);
11322d61bbb3SSatish Balay 
1133ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
1134ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
1135ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
1136ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
1137ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
1138ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
1139ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
1140ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
1141ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
1142ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
1143ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
1144ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
1145ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
1146ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
1147ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
11482d61bbb3SSatish Balay 
1149ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
1150ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
1151ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
1152ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
1153ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
1154ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
1155ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
11562d61bbb3SSatish Balay 
1157ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
1158ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
1159ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
1160ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
1161ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
1162ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_6(Mat,Vec,Vec);
1163ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
1164ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
11652d61bbb3SSatish Balay 
1166ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
1167ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
1168ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
1169ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
1170ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
1171ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
1172ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
1173ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
11742d61bbb3SSatish Balay 
11754a2ae208SSatish Balay #undef __FUNCT__
11764a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ"
11775ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
11782d61bbb3SSatish Balay {
11792d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
11802d61bbb3SSatish Balay   Mat         outA;
11812d61bbb3SSatish Balay   int         ierr;
1182667159a5SBarry Smith   PetscTruth  row_identity,col_identity;
11832d61bbb3SSatish Balay 
11842d61bbb3SSatish Balay   PetscFunctionBegin;
118529bbc08cSBarry Smith   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
1186667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1187667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1188667159a5SBarry Smith   if (!row_identity || !col_identity) {
118929bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1190667159a5SBarry Smith   }
11912d61bbb3SSatish Balay 
11922d61bbb3SSatish Balay   outA          = inA;
11932d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
11942d61bbb3SSatish Balay 
11952d61bbb3SSatish Balay   if (!a->diag) {
1196c4992f7dSBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
11972d61bbb3SSatish Balay   }
11982d61bbb3SSatish Balay   /*
119915091d37SBarry Smith       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
120015091d37SBarry Smith       for ILU(0) factorization with natural ordering
12012d61bbb3SSatish Balay   */
120215091d37SBarry Smith   switch (a->bs) {
1203f1af5d2fSBarry Smith   case 1:
1204e37c484bSBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1205e37c484bSBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_1_NaturalOrdering;
1206e37c484bSBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
1207b0a32e0cSBarry Smith     PetscLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
120815091d37SBarry Smith   case 2:
120915091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
121015091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_2_NaturalOrdering;
12117c922b88SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
1212b0a32e0cSBarry Smith     PetscLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
121315091d37SBarry Smith     break;
121415091d37SBarry Smith   case 3:
121515091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
121615091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_3_NaturalOrdering;
12177c922b88SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
1218b0a32e0cSBarry Smith     PetscLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
121915091d37SBarry Smith     break;
122015091d37SBarry Smith   case 4:
122135662bc6SKris Buschelman #if defined(PETSC_HAVE_SSE)
122235662bc6SKris Buschelman     PetscTruth sse_enabled;
122335662bc6SKris Buschelman     ierr = PetscSSEIsEnabled(&sse_enabled);CHKERRQ(ierr);
122435662bc6SKris Buschelman     if (sse_enabled) {
122535662bc6SKris Buschelman       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE;
122635662bc6SKris Buschelman     } else {
1227667159a5SBarry Smith       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
122835662bc6SKris Buschelman     }
122935662bc6SKris Buschelman #else
123035662bc6SKris Buschelman     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
123135662bc6SKris Buschelman #endif
1232f830108cSBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
12337c922b88SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering;
1234b0a32e0cSBarry Smith     PetscLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
123515091d37SBarry Smith     break;
123615091d37SBarry Smith   case 5:
1237667159a5SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
1238667159a5SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_5_NaturalOrdering;
12397c922b88SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
1240b0a32e0cSBarry Smith     PetscLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
124115091d37SBarry Smith     break;
124215091d37SBarry Smith   case 6:
124315091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
124415091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_6_NaturalOrdering;
12457c922b88SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
1246b0a32e0cSBarry Smith     PetscLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
124715091d37SBarry Smith     break;
124815091d37SBarry Smith   case 7:
124915091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
12507c922b88SBarry Smith     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering;
125115091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_7_NaturalOrdering;
1252b0a32e0cSBarry Smith     PetscLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
125315091d37SBarry Smith     break;
1254c38d4ed2SBarry Smith   default:
1255c38d4ed2SBarry Smith     a->row        = row;
1256c38d4ed2SBarry Smith     a->col        = col;
1257c38d4ed2SBarry Smith     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1258c38d4ed2SBarry Smith     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1259c38d4ed2SBarry Smith 
1260c38d4ed2SBarry Smith     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
12614c49b128SBarry Smith     ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr);
1262b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
1263c38d4ed2SBarry Smith 
1264c38d4ed2SBarry Smith     if (!a->solve_work) {
1265b0a32e0cSBarry Smith       ierr = PetscMalloc((inA->m+a->bs)*sizeof(Scalar),&a->solve_work);CHKERRQ(ierr);
1266b0a32e0cSBarry Smith       PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(Scalar));
1267c38d4ed2SBarry Smith     }
12682d61bbb3SSatish Balay   }
12692d61bbb3SSatish Balay 
1270667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1271667159a5SBarry Smith 
12722d61bbb3SSatish Balay   PetscFunctionReturn(0);
12732d61bbb3SSatish Balay }
12744a2ae208SSatish Balay #undef __FUNCT__
12754a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ"
1276ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1277ba4ca20aSSatish Balay {
1278c38d4ed2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1279ba4ca20aSSatish Balay   MPI_Comm          comm = A->comm;
1280d132466eSBarry Smith   int               ierr;
1281ba4ca20aSSatish Balay 
12823a40ed3dSBarry Smith   PetscFunctionBegin;
1283c38d4ed2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1284d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1285d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
12863a40ed3dSBarry Smith   PetscFunctionReturn(0);
1287ba4ca20aSSatish Balay }
1288d9b7c43dSSatish Balay 
1289fb2e594dSBarry Smith EXTERN_C_BEGIN
12904a2ae208SSatish Balay #undef __FUNCT__
12914a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
129227a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
129327a8da17SBarry Smith {
129427a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
129527a8da17SBarry Smith   int         i,nz,n;
129627a8da17SBarry Smith 
129727a8da17SBarry Smith   PetscFunctionBegin;
129827a8da17SBarry Smith   nz = baij->maxnz;
1299273d9f13SBarry Smith   n  = mat->n;
130027a8da17SBarry Smith   for (i=0; i<nz; i++) {
130127a8da17SBarry Smith     baij->j[i] = indices[i];
130227a8da17SBarry Smith   }
130327a8da17SBarry Smith   baij->nz = nz;
130427a8da17SBarry Smith   for (i=0; i<n; i++) {
130527a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
130627a8da17SBarry Smith   }
130727a8da17SBarry Smith 
130827a8da17SBarry Smith   PetscFunctionReturn(0);
130927a8da17SBarry Smith }
1310fb2e594dSBarry Smith EXTERN_C_END
131127a8da17SBarry Smith 
13124a2ae208SSatish Balay #undef __FUNCT__
13134a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices"
131427a8da17SBarry Smith /*@
131527a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
131627a8da17SBarry Smith        in the matrix.
131727a8da17SBarry Smith 
131827a8da17SBarry Smith   Input Parameters:
131927a8da17SBarry Smith +  mat - the SeqBAIJ matrix
132027a8da17SBarry Smith -  indices - the column indices
132127a8da17SBarry Smith 
132215091d37SBarry Smith   Level: advanced
132315091d37SBarry Smith 
132427a8da17SBarry Smith   Notes:
132527a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
132627a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
132727a8da17SBarry Smith   of the MatSetValues() operation.
132827a8da17SBarry Smith 
132927a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
133027a8da17SBarry Smith   MatCreateSeqBAIJ().
133127a8da17SBarry Smith 
133227a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
133327a8da17SBarry Smith 
133427a8da17SBarry Smith @*/
133527a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
133627a8da17SBarry Smith {
133727a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
133827a8da17SBarry Smith 
133927a8da17SBarry Smith   PetscFunctionBegin;
134027a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1341b9617806SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)())&f);CHKERRQ(ierr);
134227a8da17SBarry Smith   if (f) {
134327a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
134427a8da17SBarry Smith   } else {
134529bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
134627a8da17SBarry Smith   }
134727a8da17SBarry Smith   PetscFunctionReturn(0);
134827a8da17SBarry Smith }
134927a8da17SBarry Smith 
13504a2ae208SSatish Balay #undef __FUNCT__
13514a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ"
1352273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v)
1353273d9f13SBarry Smith {
1354273d9f13SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
1355273d9f13SBarry Smith   int          ierr,i,j,n,row,bs,*ai,*aj,mbs;
1356273d9f13SBarry Smith   PetscReal    atmp;
1357273d9f13SBarry Smith   Scalar       *x,zero = 0.0;
1358273d9f13SBarry Smith   MatScalar    *aa;
1359273d9f13SBarry Smith   int          ncols,brow,krow,kcol;
1360273d9f13SBarry Smith 
1361273d9f13SBarry Smith   PetscFunctionBegin;
1362273d9f13SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1363273d9f13SBarry Smith   bs   = a->bs;
1364273d9f13SBarry Smith   aa   = a->a;
1365273d9f13SBarry Smith   ai   = a->i;
1366273d9f13SBarry Smith   aj   = a->j;
1367273d9f13SBarry Smith   mbs = a->mbs;
1368273d9f13SBarry Smith 
1369273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1370273d9f13SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1371273d9f13SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1372273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1373273d9f13SBarry Smith   for (i=0; i<mbs; i++) {
1374273d9f13SBarry Smith     ncols = ai[1] - ai[0]; ai++;
1375273d9f13SBarry Smith     brow  = bs*i;
1376273d9f13SBarry Smith     for (j=0; j<ncols; j++){
1377273d9f13SBarry Smith       /* bcol = bs*(*aj); */
1378273d9f13SBarry Smith       for (kcol=0; kcol<bs; kcol++){
1379273d9f13SBarry Smith         for (krow=0; krow<bs; krow++){
1380273d9f13SBarry Smith           atmp = PetscAbsScalar(*aa); aa++;
1381273d9f13SBarry Smith           row = brow + krow;    /* row index */
1382273d9f13SBarry Smith           /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */
1383273d9f13SBarry Smith           if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp;
1384273d9f13SBarry Smith         }
1385273d9f13SBarry Smith       }
1386273d9f13SBarry Smith       aj++;
1387273d9f13SBarry Smith     }
1388273d9f13SBarry Smith   }
1389273d9f13SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1390273d9f13SBarry Smith   PetscFunctionReturn(0);
1391273d9f13SBarry Smith }
1392273d9f13SBarry Smith 
13934a2ae208SSatish Balay #undef __FUNCT__
13944a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ"
1395273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A)
1396273d9f13SBarry Smith {
1397273d9f13SBarry Smith   int        ierr;
1398273d9f13SBarry Smith 
1399273d9f13SBarry Smith   PetscFunctionBegin;
1400273d9f13SBarry Smith   ierr =  MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1401273d9f13SBarry Smith   PetscFunctionReturn(0);
1402273d9f13SBarry Smith }
1403273d9f13SBarry Smith 
14044a2ae208SSatish Balay #undef __FUNCT__
14054a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ"
1406f2a5309cSSatish Balay int MatGetArray_SeqBAIJ(Mat A,Scalar **array)
1407f2a5309cSSatish Balay {
1408f2a5309cSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1409f2a5309cSSatish Balay   PetscFunctionBegin;
1410f2a5309cSSatish Balay   *array = a->a;
1411f2a5309cSSatish Balay   PetscFunctionReturn(0);
1412f2a5309cSSatish Balay }
1413f2a5309cSSatish Balay 
14144a2ae208SSatish Balay #undef __FUNCT__
14154a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ"
1416f2a5309cSSatish Balay int MatRestoreArray_SeqBAIJ(Mat A,Scalar **array)
1417f2a5309cSSatish Balay {
1418f2a5309cSSatish Balay   PetscFunctionBegin;
1419f2a5309cSSatish Balay   PetscFunctionReturn(0);
1420f2a5309cSSatish Balay }
1421f2a5309cSSatish Balay 
14222593348eSBarry Smith /* -------------------------------------------------------------------*/
1423cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1424cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1425cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1426cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1427cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
14287c922b88SBarry Smith        MatMultTranspose_SeqBAIJ,
14297c922b88SBarry Smith        MatMultTransposeAdd_SeqBAIJ,
1430cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1431cc2dc46cSBarry Smith        0,
1432cc2dc46cSBarry Smith        0,
1433cc2dc46cSBarry Smith        0,
1434cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1435cc2dc46cSBarry Smith        0,
1436b6490206SBarry Smith        0,
1437f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1438cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1439cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1440cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1441cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1442cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1443b6490206SBarry Smith        0,
1444cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1445cc2dc46cSBarry Smith        0,
1446cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1447cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1448cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1449cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1450cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1451cc2dc46cSBarry Smith        0,
1452cc2dc46cSBarry Smith        0,
1453273d9f13SBarry Smith        MatSetUpPreallocation_SeqBAIJ,
1454273d9f13SBarry Smith        0,
1455cc2dc46cSBarry Smith        MatGetOwnershipRange_SeqBAIJ,
1456cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1457cc2dc46cSBarry Smith        0,
1458f2a5309cSSatish Balay        MatGetArray_SeqBAIJ,
1459f2a5309cSSatish Balay        MatRestoreArray_SeqBAIJ,
14602e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1461cc2dc46cSBarry Smith        0,
1462cc2dc46cSBarry Smith        0,
1463cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1464cc2dc46cSBarry Smith        0,
1465cc2dc46cSBarry Smith        0,
1466cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1467cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1468cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1469cc2dc46cSBarry Smith        0,
1470cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1471cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1472cc2dc46cSBarry Smith        0,
1473cc2dc46cSBarry Smith        0,
1474cc2dc46cSBarry Smith        0,
1475cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
14763b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
147792c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1478cc2dc46cSBarry Smith        0,
1479cc2dc46cSBarry Smith        0,
1480cc2dc46cSBarry Smith        0,
1481cc2dc46cSBarry Smith        0,
1482cc2dc46cSBarry Smith        0,
1483cc2dc46cSBarry Smith        0,
1484d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1485cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1486b9b97703SBarry Smith        MatDestroy_SeqBAIJ,
1487b9b97703SBarry Smith        MatView_SeqBAIJ,
1488273d9f13SBarry Smith        MatGetMaps_Petsc,
1489273d9f13SBarry Smith        0,
1490273d9f13SBarry Smith        0,
1491273d9f13SBarry Smith        0,
1492273d9f13SBarry Smith        0,
1493273d9f13SBarry Smith        0,
1494273d9f13SBarry Smith        0,
1495273d9f13SBarry Smith        MatGetRowMax_SeqBAIJ,
1496273d9f13SBarry Smith        MatConvert_Basic};
14972593348eSBarry Smith 
14983e90b805SBarry Smith EXTERN_C_BEGIN
14994a2ae208SSatish Balay #undef __FUNCT__
15004a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ"
15013e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
15023e90b805SBarry Smith {
15033e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1504273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1505d132466eSBarry Smith   int         ierr;
15063e90b805SBarry Smith 
15073e90b805SBarry Smith   PetscFunctionBegin;
15083e90b805SBarry Smith   if (aij->nonew != 1) {
150929bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15103e90b805SBarry Smith   }
15113e90b805SBarry Smith 
15123e90b805SBarry Smith   /* allocate space for values if not already there */
15133e90b805SBarry Smith   if (!aij->saved_values) {
1514f2a5309cSSatish Balay     ierr = PetscMalloc((nz+1)*sizeof(Scalar),&aij->saved_values);CHKERRQ(ierr);
15153e90b805SBarry Smith   }
15163e90b805SBarry Smith 
15173e90b805SBarry Smith   /* copy values over */
1518d132466eSBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr);
15193e90b805SBarry Smith   PetscFunctionReturn(0);
15203e90b805SBarry Smith }
15213e90b805SBarry Smith EXTERN_C_END
15223e90b805SBarry Smith 
15233e90b805SBarry Smith EXTERN_C_BEGIN
15244a2ae208SSatish Balay #undef __FUNCT__
15254a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ"
152632e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
15273e90b805SBarry Smith {
15283e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1529273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
15303e90b805SBarry Smith 
15313e90b805SBarry Smith   PetscFunctionBegin;
15323e90b805SBarry Smith   if (aij->nonew != 1) {
153329bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15343e90b805SBarry Smith   }
15353e90b805SBarry Smith   if (!aij->saved_values) {
153629bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
15373e90b805SBarry Smith   }
15383e90b805SBarry Smith 
15393e90b805SBarry Smith   /* copy values over */
1540549d3d68SSatish Balay   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr);
15413e90b805SBarry Smith   PetscFunctionReturn(0);
15423e90b805SBarry Smith }
15433e90b805SBarry Smith EXTERN_C_END
15443e90b805SBarry Smith 
1545273d9f13SBarry Smith EXTERN_C_BEGIN
1546273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*);
1547273d9f13SBarry Smith EXTERN_C_END
1548273d9f13SBarry Smith 
1549273d9f13SBarry Smith EXTERN_C_BEGIN
15504a2ae208SSatish Balay #undef __FUNCT__
15514a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ"
1552273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B)
15532593348eSBarry Smith {
1554273d9f13SBarry Smith   int         ierr,size;
1555b6490206SBarry Smith   Mat_SeqBAIJ *b;
15563b2fbd54SBarry Smith 
15573a40ed3dSBarry Smith   PetscFunctionBegin;
1558273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
155929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1560b6490206SBarry Smith 
1561273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1562273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1563b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr);
1564b0a32e0cSBarry Smith   B->data = (void*)b;
1565549d3d68SSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1566549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
15672593348eSBarry Smith   B->factor           = 0;
15682593348eSBarry Smith   B->lupivotthreshold = 1.0;
156990f02eecSBarry Smith   B->mapping          = 0;
15702593348eSBarry Smith   b->row              = 0;
15712593348eSBarry Smith   b->col              = 0;
1572e51c0b9cSSatish Balay   b->icol             = 0;
15732593348eSBarry Smith   b->reallocs         = 0;
15743e90b805SBarry Smith   b->saved_values     = 0;
1575cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
1576563b5814SBarry Smith   b->setvalueslen     = 0;
1577563b5814SBarry Smith   b->setvaluescopy    = PETSC_NULL;
1578563b5814SBarry Smith #endif
15792593348eSBarry Smith 
1580273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
1581273d9f13SBarry Smith   ierr = MapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1582a5ae1ecdSBarry Smith 
1583c4992f7dSBarry Smith   b->sorted           = PETSC_FALSE;
1584c4992f7dSBarry Smith   b->roworiented      = PETSC_TRUE;
15852593348eSBarry Smith   b->nonew            = 0;
15862593348eSBarry Smith   b->diag             = 0;
15872593348eSBarry Smith   b->solve_work       = 0;
1588de6a44a3SBarry Smith   b->mult_work        = 0;
15892593348eSBarry Smith   b->spptr            = 0;
15900e6d2581SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1591883fce79SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
15924e220ebcSLois Curfman McInnes 
1593f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
15943e90b805SBarry Smith                                      "MatStoreValues_SeqBAIJ",
1595bc4b532fSSatish Balay                                       MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1596f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
15973e90b805SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
1598bc4b532fSSatish Balay                                       MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1599f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
160027a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
1601bc4b532fSSatish Balay                                       MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
1602273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_SeqBAIJ_SeqAIJ_C",
1603273d9f13SBarry Smith                                      "MatConvert_SeqBAIJ_SeqAIJ",
1604273d9f13SBarry Smith                                       MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr);
16053a40ed3dSBarry Smith   PetscFunctionReturn(0);
16062593348eSBarry Smith }
1607273d9f13SBarry Smith EXTERN_C_END
16082593348eSBarry Smith 
16094a2ae208SSatish Balay #undef __FUNCT__
16104a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ"
16112e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
16122593348eSBarry Smith {
16132593348eSBarry Smith   Mat         C;
1614b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data;
161527a8da17SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1616de6a44a3SBarry Smith 
16173a40ed3dSBarry Smith   PetscFunctionBegin;
161829bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
16192593348eSBarry Smith 
16202593348eSBarry Smith   *B = 0;
1621273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1622273d9f13SBarry Smith   ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr);
1623273d9f13SBarry Smith   c    = (Mat_SeqBAIJ*)C->data;
162444cd7ae7SLois Curfman McInnes 
1625b6490206SBarry Smith   c->bs         = a->bs;
16269df24120SSatish Balay   c->bs2        = a->bs2;
1627b6490206SBarry Smith   c->mbs        = a->mbs;
1628b6490206SBarry Smith   c->nbs        = a->nbs;
162935d8aa7fSBarry Smith   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
16302593348eSBarry Smith 
1631b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1632b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
1633b6490206SBarry Smith   for (i=0; i<mbs; i++) {
16342593348eSBarry Smith     c->imax[i] = a->imax[i];
16352593348eSBarry Smith     c->ilen[i] = a->ilen[i];
16362593348eSBarry Smith   }
16372593348eSBarry Smith 
16382593348eSBarry Smith   /* allocate the matrix space */
1639c4992f7dSBarry Smith   c->singlemalloc = PETSC_TRUE;
16403f1db9ecSBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1641b0a32e0cSBarry Smith   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
16427e67e3f9SSatish Balay   c->j = (int*)(c->a + nz*bs2);
1643de6a44a3SBarry Smith   c->i = c->j + nz;
1644549d3d68SSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1645b6490206SBarry Smith   if (mbs > 0) {
1646549d3d68SSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
16472e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1648549d3d68SSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
16492e8a6d31SBarry Smith     } else {
1650549d3d68SSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
16512593348eSBarry Smith     }
16522593348eSBarry Smith   }
16532593348eSBarry Smith 
1654b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
16552593348eSBarry Smith   c->sorted      = a->sorted;
16562593348eSBarry Smith   c->roworiented = a->roworiented;
16572593348eSBarry Smith   c->nonew       = a->nonew;
16582593348eSBarry Smith 
16592593348eSBarry Smith   if (a->diag) {
1660b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1661b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
1662b6490206SBarry Smith     for (i=0; i<mbs; i++) {
16632593348eSBarry Smith       c->diag[i] = a->diag[i];
16642593348eSBarry Smith     }
166598305bb5SBarry Smith   } else c->diag        = 0;
16662593348eSBarry Smith   c->nz                 = a->nz;
16672593348eSBarry Smith   c->maxnz              = a->maxnz;
16682593348eSBarry Smith   c->solve_work         = 0;
16692593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
16707fc0212eSBarry Smith   c->mult_work          = 0;
1671273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
1672273d9f13SBarry Smith   C->assembled          = PETSC_TRUE;
16732593348eSBarry Smith   *B = C;
1674b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
16753a40ed3dSBarry Smith   PetscFunctionReturn(0);
16762593348eSBarry Smith }
16772593348eSBarry Smith 
1678273d9f13SBarry Smith EXTERN_C_BEGIN
16794a2ae208SSatish Balay #undef __FUNCT__
16804a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ"
1681b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A)
16822593348eSBarry Smith {
1683b6490206SBarry Smith   Mat_SeqBAIJ  *a;
16842593348eSBarry Smith   Mat          B;
1685f1af5d2fSBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1686b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
168735aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1688a7c10996SSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
1689b6490206SBarry Smith   Scalar       *aa;
169019bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
16912593348eSBarry Smith 
16923a40ed3dSBarry Smith   PetscFunctionBegin;
1693b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1694de6a44a3SBarry Smith   bs2  = bs*bs;
1695b6490206SBarry Smith 
1696d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
169729bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1698b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
16990752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
170029bbc08cSBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
17012593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
17022593348eSBarry Smith 
1703d64ed03dSBarry Smith   if (header[3] < 0) {
170429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ");
1705d64ed03dSBarry Smith   }
1706d64ed03dSBarry Smith 
170729bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
170835aab85fSBarry Smith 
170935aab85fSBarry Smith   /*
171035aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
171135aab85fSBarry Smith     divisible by the blocksize
171235aab85fSBarry Smith   */
1713b6490206SBarry Smith   mbs        = M/bs;
171435aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
171535aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
171635aab85fSBarry Smith   else                  mbs++;
171735aab85fSBarry Smith   if (extra_rows) {
1718b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
171935aab85fSBarry Smith   }
1720b6490206SBarry Smith 
17212593348eSBarry Smith   /* read in row lengths */
1722b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
17230752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
172435aab85fSBarry Smith   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
17252593348eSBarry Smith 
1726b6490206SBarry Smith   /* read in column indices */
1727b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
17280752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
172935aab85fSBarry Smith   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
1730b6490206SBarry Smith 
1731b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1732b0a32e0cSBarry Smith   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
1733549d3d68SSatish Balay   ierr     = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
1734b0a32e0cSBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
1735549d3d68SSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
173635aab85fSBarry Smith   masked   = mask + mbs;
1737b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1738b6490206SBarry Smith   for (i=0; i<mbs; i++) {
173935aab85fSBarry Smith     nmask = 0;
1740b6490206SBarry Smith     for (j=0; j<bs; j++) {
1741b6490206SBarry Smith       kmax = rowlengths[rowcount];
1742b6490206SBarry Smith       for (k=0; k<kmax; k++) {
174335aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
174435aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1745b6490206SBarry Smith       }
1746b6490206SBarry Smith       rowcount++;
1747b6490206SBarry Smith     }
174835aab85fSBarry Smith     browlengths[i] += nmask;
174935aab85fSBarry Smith     /* zero out the mask elements we set */
175035aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1751b6490206SBarry Smith   }
1752b6490206SBarry Smith 
17532593348eSBarry Smith   /* create our matrix */
17543eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
17552593348eSBarry Smith   B = *A;
1756b6490206SBarry Smith   a = (Mat_SeqBAIJ*)B->data;
17572593348eSBarry Smith 
17582593348eSBarry Smith   /* set matrix "i" values */
1759de6a44a3SBarry Smith   a->i[0] = 0;
1760b6490206SBarry Smith   for (i=1; i<= mbs; i++) {
1761b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1762b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
17632593348eSBarry Smith   }
1764b6490206SBarry Smith   a->nz         = 0;
1765b6490206SBarry Smith   for (i=0; i<mbs; i++) a->nz += browlengths[i];
17662593348eSBarry Smith 
1767b6490206SBarry Smith   /* read in nonzero values */
1768b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(Scalar),&aa);CHKERRQ(ierr);
17690752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
177035aab85fSBarry Smith   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
1771b6490206SBarry Smith 
1772b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1773b6490206SBarry Smith   nzcount = 0; jcount = 0;
1774b6490206SBarry Smith   for (i=0; i<mbs; i++) {
1775b6490206SBarry Smith     nzcountb = nzcount;
177635aab85fSBarry Smith     nmask    = 0;
1777b6490206SBarry Smith     for (j=0; j<bs; j++) {
1778b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1779b6490206SBarry Smith       for (k=0; k<kmax; k++) {
178035aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
178135aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1782b6490206SBarry Smith       }
1783b6490206SBarry Smith     }
1784de6a44a3SBarry Smith     /* sort the masked values */
1785433994e6SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1786de6a44a3SBarry Smith 
1787b6490206SBarry Smith     /* set "j" values into matrix */
1788b6490206SBarry Smith     maskcount = 1;
178935aab85fSBarry Smith     for (j=0; j<nmask; j++) {
179035aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1791de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1792b6490206SBarry Smith     }
1793b6490206SBarry Smith     /* set "a" values into matrix */
1794de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1795b6490206SBarry Smith     for (j=0; j<bs; j++) {
1796b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1797b6490206SBarry Smith       for (k=0; k<kmax; k++) {
1798de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1799de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1800de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1801de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1802375fe846SBarry Smith         a->a[idx] = (MatScalar)aa[nzcountb++];
1803b6490206SBarry Smith       }
1804b6490206SBarry Smith     }
180535aab85fSBarry Smith     /* zero out the mask elements we set */
180635aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1807b6490206SBarry Smith   }
180829bbc08cSBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
1809b6490206SBarry Smith 
1810606d414cSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1811606d414cSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1812606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
1813606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
1814606d414cSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
1815b6490206SBarry Smith 
1816b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1817de6a44a3SBarry Smith 
18189c01be13SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
18193a40ed3dSBarry Smith   PetscFunctionReturn(0);
18202593348eSBarry Smith }
1821273d9f13SBarry Smith EXTERN_C_END
18222593348eSBarry Smith 
18234a2ae208SSatish Balay #undef __FUNCT__
18244a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ"
1825273d9f13SBarry Smith /*@C
1826273d9f13SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
1827273d9f13SBarry Smith    compressed row) format.  For good matrix assembly performance the
1828273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1829273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1830273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
18312593348eSBarry Smith 
1832273d9f13SBarry Smith    Collective on MPI_Comm
1833273d9f13SBarry Smith 
1834273d9f13SBarry Smith    Input Parameters:
1835273d9f13SBarry Smith +  comm - MPI communicator, set to PETSC_COMM_SELF
1836273d9f13SBarry Smith .  bs - size of block
1837273d9f13SBarry Smith .  m - number of rows
1838273d9f13SBarry Smith .  n - number of columns
183935d8aa7fSBarry Smith .  nz - number of nonzero blocks  per block row (same for all rows)
184035d8aa7fSBarry Smith -  nnz - array containing the number of nonzero blocks in the various block rows
1841273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1842273d9f13SBarry Smith 
1843273d9f13SBarry Smith    Output Parameter:
1844273d9f13SBarry Smith .  A - the matrix
1845273d9f13SBarry Smith 
1846273d9f13SBarry Smith    Options Database Keys:
1847273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1848273d9f13SBarry Smith                      block calculations (much slower)
1849273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1850273d9f13SBarry Smith 
1851273d9f13SBarry Smith    Level: intermediate
1852273d9f13SBarry Smith 
1853273d9f13SBarry Smith    Notes:
185435d8aa7fSBarry Smith    A nonzero block is any block that as 1 or more nonzeros in it
185535d8aa7fSBarry Smith 
1856273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1857273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1858273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1859273d9f13SBarry Smith 
1860273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1861273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1862273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1863273d9f13SBarry Smith    matrices.
1864273d9f13SBarry Smith 
1865273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1866273d9f13SBarry Smith @*/
1867273d9f13SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1868273d9f13SBarry Smith {
1869273d9f13SBarry Smith   int ierr;
1870273d9f13SBarry Smith 
1871273d9f13SBarry Smith   PetscFunctionBegin;
1872273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1873273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr);
1874273d9f13SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1875273d9f13SBarry Smith   PetscFunctionReturn(0);
1876273d9f13SBarry Smith }
1877273d9f13SBarry Smith 
18784a2ae208SSatish Balay #undef __FUNCT__
18794a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation"
1880273d9f13SBarry Smith /*@C
1881273d9f13SBarry Smith    MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros
1882273d9f13SBarry Smith    per row in the matrix. For good matrix assembly performance the
1883273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1884273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1885273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
1886273d9f13SBarry Smith 
1887273d9f13SBarry Smith    Collective on MPI_Comm
1888273d9f13SBarry Smith 
1889273d9f13SBarry Smith    Input Parameters:
1890273d9f13SBarry Smith +  A - the matrix
1891273d9f13SBarry Smith .  bs - size of block
1892273d9f13SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1893273d9f13SBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
1894273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1895273d9f13SBarry Smith 
1896273d9f13SBarry Smith    Options Database Keys:
1897273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1898273d9f13SBarry Smith                      block calculations (much slower)
1899273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1900273d9f13SBarry Smith 
1901273d9f13SBarry Smith    Level: intermediate
1902273d9f13SBarry Smith 
1903273d9f13SBarry Smith    Notes:
1904273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1905273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1906273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1907273d9f13SBarry Smith 
1908273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1909273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1910273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1911273d9f13SBarry Smith    matrices.
1912273d9f13SBarry Smith 
1913273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1914273d9f13SBarry Smith @*/
1915273d9f13SBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
1916273d9f13SBarry Smith {
1917273d9f13SBarry Smith   Mat_SeqBAIJ *b;
1918273d9f13SBarry Smith   int         i,len,ierr,mbs,nbs,bs2,newbs = bs;
1919273d9f13SBarry Smith   PetscTruth  flg;
1920273d9f13SBarry Smith 
1921273d9f13SBarry Smith   PetscFunctionBegin;
1922273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flg);CHKERRQ(ierr);
1923273d9f13SBarry Smith   if (!flg) PetscFunctionReturn(0);
1924273d9f13SBarry Smith 
1925273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1926b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr);
1927273d9f13SBarry Smith   if (nnz && newbs != bs) {
1928273d9f13SBarry Smith     SETERRQ(1,"Cannot change blocksize from command line if setting nnz");
1929273d9f13SBarry Smith   }
1930273d9f13SBarry Smith   bs   = newbs;
1931273d9f13SBarry Smith 
1932273d9f13SBarry Smith   mbs  = B->m/bs;
1933273d9f13SBarry Smith   nbs  = B->n/bs;
1934273d9f13SBarry Smith   bs2  = bs*bs;
1935273d9f13SBarry Smith 
1936273d9f13SBarry Smith   if (mbs*bs!=B->m || nbs*bs!=B->n) {
193735d8aa7fSBarry Smith     SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs);
1938273d9f13SBarry Smith   }
1939273d9f13SBarry Smith 
1940435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1941435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
1942273d9f13SBarry Smith   if (nnz) {
1943273d9f13SBarry Smith     for (i=0; i<mbs; i++) {
1944273d9f13SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
19453a7fca6bSBarry Smith       if (nnz[i] > nbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %d value %d rowlength %d",i,nnz[i],nbs);
1946273d9f13SBarry Smith     }
1947273d9f13SBarry Smith   }
1948273d9f13SBarry Smith 
1949273d9f13SBarry Smith   b       = (Mat_SeqBAIJ*)B->data;
1950b0a32e0cSBarry Smith   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
1951273d9f13SBarry Smith   if (!flg) {
1952273d9f13SBarry Smith     switch (bs) {
1953273d9f13SBarry Smith     case 1:
1954273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1955273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1956273d9f13SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_1;
1957273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1958273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
1959273d9f13SBarry Smith       break;
1960273d9f13SBarry Smith     case 2:
1961273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1962273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1963273d9f13SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_2;
1964273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1965273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
1966273d9f13SBarry Smith       break;
1967273d9f13SBarry Smith     case 3:
1968273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1969273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1970273d9f13SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_3;
1971273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1972273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
1973273d9f13SBarry Smith       break;
1974273d9f13SBarry Smith     case 4:
1975273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1976273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1977273d9f13SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_4;
1978273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1979273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
1980273d9f13SBarry Smith       break;
1981273d9f13SBarry Smith     case 5:
1982273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1983273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1984273d9f13SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_5;
1985273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1986273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
1987273d9f13SBarry Smith       break;
1988273d9f13SBarry Smith     case 6:
1989273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
1990273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_6;
1991273d9f13SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_6;
1992273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
1993273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
1994273d9f13SBarry Smith       break;
1995273d9f13SBarry Smith     case 7:
1996273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1997273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1998273d9f13SBarry Smith       B->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_7;
1999273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
2000273d9f13SBarry Smith       break;
2001273d9f13SBarry Smith     default:
2002273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_N;
2003273d9f13SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_N;
2004273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_N;
2005273d9f13SBarry Smith       break;
2006273d9f13SBarry Smith     }
2007273d9f13SBarry Smith   }
2008273d9f13SBarry Smith   b->bs      = bs;
2009273d9f13SBarry Smith   b->mbs     = mbs;
2010273d9f13SBarry Smith   b->nbs     = nbs;
2011b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
2012273d9f13SBarry Smith   if (!nnz) {
2013435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2014273d9f13SBarry Smith     else if (nz <= 0)        nz = 1;
2015273d9f13SBarry Smith     for (i=0; i<mbs; i++) b->imax[i] = nz;
2016273d9f13SBarry Smith     nz = nz*mbs;
2017273d9f13SBarry Smith   } else {
2018273d9f13SBarry Smith     nz = 0;
2019273d9f13SBarry Smith     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2020273d9f13SBarry Smith   }
2021273d9f13SBarry Smith 
2022273d9f13SBarry Smith   /* allocate the matrix space */
2023273d9f13SBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
2024b0a32e0cSBarry Smith   ierr  = PetscMalloc(len,&b->a);CHKERRQ(ierr);
2025273d9f13SBarry Smith   ierr  = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
2026273d9f13SBarry Smith   b->j  = (int*)(b->a + nz*bs2);
2027273d9f13SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
2028273d9f13SBarry Smith   b->i  = b->j + nz;
2029273d9f13SBarry Smith   b->singlemalloc = PETSC_TRUE;
2030273d9f13SBarry Smith 
2031273d9f13SBarry Smith   b->i[0] = 0;
2032273d9f13SBarry Smith   for (i=1; i<mbs+1; i++) {
2033273d9f13SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
2034273d9f13SBarry Smith   }
2035273d9f13SBarry Smith 
2036273d9f13SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
203716d6aa21SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
2038b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
2039273d9f13SBarry Smith   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
2040273d9f13SBarry Smith 
2041273d9f13SBarry Smith   b->bs               = bs;
2042273d9f13SBarry Smith   b->bs2              = bs2;
2043273d9f13SBarry Smith   b->mbs              = mbs;
2044273d9f13SBarry Smith   b->nz               = 0;
2045273d9f13SBarry Smith   b->maxnz            = nz*bs2;
2046273d9f13SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
2047273d9f13SBarry Smith   PetscFunctionReturn(0);
2048273d9f13SBarry Smith }
20492593348eSBarry Smith 
2050