xref: /petsc/src/mat/impls/baij/seq/baij.c (revision a56a16c87eec4e97c4cf0ce1e098241a17b4c30e)
1bba1ac68SSatish Balay /*$Id: baij.c,v 1.245 2001/08/07 03:02:55 balay Exp bsmith $*/
22593348eSBarry Smith 
32593348eSBarry Smith /*
4b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
52593348eSBarry Smith   matrix storage format.
62593348eSBarry Smith */
770f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
81a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
91a6a6d98SBarry Smith #include "src/inline/spops.h"
10325e03aeSBarry Smith #include "petscsys.h"                     /*I "petscmat.h" I*/
113b547af2SSatish Balay 
1295d5f7c2SBarry Smith /*  UGLY, ugly, ugly
1387828ca2SBarry Smith    When MatScalar == PetscScalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does
143477af42SKris Buschelman    not exist. Otherwise ..._MatScalar() takes matrix dlements 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");
213bd648c37SKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
214bd648c37SKris Buschelman     if (a->bs==4) {
21594ee7fc8SKris Buschelman       PetscTruth sse_enabled_local,sse_enabled_global;
216bd648c37SKris Buschelman       int        ierr;
21794ee7fc8SKris Buschelman 
21894ee7fc8SKris Buschelman       sse_enabled_local  = PETSC_FALSE;
21994ee7fc8SKris Buschelman       sse_enabled_global = PETSC_FALSE;
22094ee7fc8SKris Buschelman 
22194ee7fc8SKris Buschelman       ierr = PetscSSEIsEnabled(A->comm,&sse_enabled_local,&sse_enabled_global);CHKERRQ(ierr);
222e64df4b9SKris Buschelman #if defined(PETSC_HAVE_SSE)
22394ee7fc8SKris Buschelman       if (sse_enabled_local) {
22454138f6bSKris Buschelman           a->single_precision_solves = PETSC_TRUE;
22554138f6bSKris Buschelman           A->ops->solve              = MatSolve_SeqBAIJ_Update;
22654138f6bSKris Buschelman           A->ops->solvetranspose     = MatSolveTranspose_SeqBAIJ_Update;
227cf242676SKris Buschelman           PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES set\n");
22854138f6bSKris Buschelman           break;
22994ee7fc8SKris Buschelman       } else {
23094ee7fc8SKris Buschelman         PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n");
2318661488fSKris Buschelman       }
232e64df4b9SKris Buschelman #else
233bd648c37SKris Buschelman       PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n");
234e64df4b9SKris Buschelman #endif
235bd648c37SKris Buschelman     }
236bd648c37SKris Buschelman     break;
237aa275fccSKris Buschelman   default:
23829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
2392d61bbb3SSatish Balay   }
2402d61bbb3SSatish Balay   PetscFunctionReturn(0);
2412d61bbb3SSatish Balay }
2422d61bbb3SSatish Balay 
2434a2ae208SSatish Balay #undef __FUNCT__
2444a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ"
24587828ca2SBarry Smith int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
2462d61bbb3SSatish Balay {
2472d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
24882502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr;
2493f1db9ecSBarry Smith   MatScalar    *aa,*aa_i;
25087828ca2SBarry Smith   PetscScalar  *v_i;
2512d61bbb3SSatish Balay 
2522d61bbb3SSatish Balay   PetscFunctionBegin;
2532d61bbb3SSatish Balay   bs  = a->bs;
2542d61bbb3SSatish Balay   ai  = a->i;
2552d61bbb3SSatish Balay   aj  = a->j;
2562d61bbb3SSatish Balay   aa  = a->a;
2572d61bbb3SSatish Balay   bs2 = a->bs2;
2582d61bbb3SSatish Balay 
259273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
2602d61bbb3SSatish Balay 
2612d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2622d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2632d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2642d61bbb3SSatish Balay   *nz = bs*M;
2652d61bbb3SSatish Balay 
2662d61bbb3SSatish Balay   if (v) {
2672d61bbb3SSatish Balay     *v = 0;
2682d61bbb3SSatish Balay     if (*nz) {
26987828ca2SBarry Smith       ierr = PetscMalloc((*nz)*sizeof(PetscScalar),v);CHKERRQ(ierr);
2702d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
2712d61bbb3SSatish Balay         v_i  = *v + i*bs;
2722d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2732d61bbb3SSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
2742d61bbb3SSatish Balay       }
2752d61bbb3SSatish Balay     }
2762d61bbb3SSatish Balay   }
2772d61bbb3SSatish Balay 
2782d61bbb3SSatish Balay   if (idx) {
2792d61bbb3SSatish Balay     *idx = 0;
2802d61bbb3SSatish Balay     if (*nz) {
281b0a32e0cSBarry Smith       ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr);
2822d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
2832d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2842d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2852d61bbb3SSatish Balay         for (j=0; j<bs; j++) {idx_i[j] = itmp++;}
2862d61bbb3SSatish Balay       }
2872d61bbb3SSatish Balay     }
2882d61bbb3SSatish Balay   }
2892d61bbb3SSatish Balay   PetscFunctionReturn(0);
2902d61bbb3SSatish Balay }
2912d61bbb3SSatish Balay 
2924a2ae208SSatish Balay #undef __FUNCT__
2934a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ"
29487828ca2SBarry Smith int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
2952d61bbb3SSatish Balay {
296606d414cSSatish Balay   int ierr;
297606d414cSSatish Balay 
2982d61bbb3SSatish Balay   PetscFunctionBegin;
299606d414cSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
300606d414cSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
3012d61bbb3SSatish Balay   PetscFunctionReturn(0);
3022d61bbb3SSatish Balay }
3032d61bbb3SSatish Balay 
3044a2ae208SSatish Balay #undef __FUNCT__
3054a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ"
3062d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
3072d61bbb3SSatish Balay {
3082d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
3092d61bbb3SSatish Balay   Mat         C;
3102d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
311273d9f13SBarry Smith   int         *rows,*cols,bs2=a->bs2;
31287828ca2SBarry Smith   PetscScalar *array;
3132d61bbb3SSatish Balay 
3142d61bbb3SSatish Balay   PetscFunctionBegin;
31529bbc08cSBarry Smith   if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place");
316b0a32e0cSBarry Smith   ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr);
317549d3d68SSatish Balay   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
3182d61bbb3SSatish Balay 
319375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
32087828ca2SBarry Smith   ierr = PetscMalloc(a->bs2*a->nz*sizeof(PetscScalar),&array);CHKERRQ(ierr);
32187828ca2SBarry Smith   for (i=0; i<a->bs2*a->nz; i++) array[i] = (PetscScalar)a->a[i];
322375fe846SBarry Smith #else
3233eda8832SBarry Smith   array = a->a;
324375fe846SBarry Smith #endif
325375fe846SBarry Smith 
3262d61bbb3SSatish Balay   for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1;
327273d9f13SBarry Smith   ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
328606d414cSSatish Balay   ierr = PetscFree(col);CHKERRQ(ierr);
329b0a32e0cSBarry Smith   ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr);
3302d61bbb3SSatish Balay   cols = rows + bs;
3312d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3322d61bbb3SSatish Balay     cols[0] = i*bs;
3332d61bbb3SSatish Balay     for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1;
3342d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
3352d61bbb3SSatish Balay     for (j=0; j<len; j++) {
3362d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
3372d61bbb3SSatish Balay       for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1;
3382d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
3392d61bbb3SSatish Balay       array += bs2;
3402d61bbb3SSatish Balay     }
3412d61bbb3SSatish Balay   }
342606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
343375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
344375fe846SBarry Smith   ierr = PetscFree(array);
345375fe846SBarry Smith #endif
3462d61bbb3SSatish Balay 
3472d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3482d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3492d61bbb3SSatish Balay 
350c4992f7dSBarry Smith   if (B) {
3512d61bbb3SSatish Balay     *B = C;
3522d61bbb3SSatish Balay   } else {
353273d9f13SBarry Smith     ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
3542d61bbb3SSatish Balay   }
3552d61bbb3SSatish Balay   PetscFunctionReturn(0);
3562d61bbb3SSatish Balay }
3572d61bbb3SSatish Balay 
3584a2ae208SSatish Balay #undef __FUNCT__
3594a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary"
360b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer)
3612593348eSBarry Smith {
362b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
3639df24120SSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
36487828ca2SBarry Smith   PetscScalar *aa;
365ce6f0cecSBarry Smith   FILE        *file;
3662593348eSBarry Smith 
3673a40ed3dSBarry Smith   PetscFunctionBegin;
368b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
369b0a32e0cSBarry Smith   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
370552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
3713b2fbd54SBarry Smith 
372273d9f13SBarry Smith   col_lens[1] = A->m;
373273d9f13SBarry Smith   col_lens[2] = A->n;
3747e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3752593348eSBarry Smith 
3762593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
377b6490206SBarry Smith   count = 0;
378b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
379b6490206SBarry Smith     for (j=0; j<bs; j++) {
380b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
381b6490206SBarry Smith     }
3822593348eSBarry Smith   }
383273d9f13SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
384606d414cSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
3852593348eSBarry Smith 
3862593348eSBarry Smith   /* store column indices (zero start index) */
387b0a32e0cSBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
388b6490206SBarry Smith   count = 0;
389b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
390b6490206SBarry Smith     for (j=0; j<bs; j++) {
391b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
392b6490206SBarry Smith         for (l=0; l<bs; l++) {
393b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3942593348eSBarry Smith         }
3952593348eSBarry Smith       }
396b6490206SBarry Smith     }
397b6490206SBarry Smith   }
3980752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
399606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
4002593348eSBarry Smith 
4012593348eSBarry Smith   /* store nonzero values */
40287828ca2SBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
403b6490206SBarry Smith   count = 0;
404b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
405b6490206SBarry Smith     for (j=0; j<bs; j++) {
406b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
407b6490206SBarry Smith         for (l=0; l<bs; l++) {
4087e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
409b6490206SBarry Smith         }
410b6490206SBarry Smith       }
411b6490206SBarry Smith     }
412b6490206SBarry Smith   }
4130752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
414606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
415ce6f0cecSBarry Smith 
416b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
417ce6f0cecSBarry Smith   if (file) {
418ce6f0cecSBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
419ce6f0cecSBarry Smith   }
4203a40ed3dSBarry Smith   PetscFunctionReturn(0);
4212593348eSBarry Smith }
4222593348eSBarry Smith 
4234a2ae208SSatish Balay #undef __FUNCT__
4244a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII"
425b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer)
4262593348eSBarry Smith {
427b6490206SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
428fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
429f3ef73ceSBarry Smith   PetscViewerFormat format;
4302593348eSBarry Smith 
4313a40ed3dSBarry Smith   PetscFunctionBegin;
432b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
433fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
434b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
435fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
43629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
437fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
438b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
43944cd7ae7SLois Curfman McInnes     for (i=0; i<a->mbs; i++) {
44044cd7ae7SLois Curfman McInnes       for (j=0; j<bs; j++) {
441b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
44244cd7ae7SLois Curfman McInnes         for (k=a->i[i]; k<a->i[i+1]; k++) {
44344cd7ae7SLois Curfman McInnes           for (l=0; l<bs; l++) {
444aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4450e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
446b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %gi",bs*a->j[k]+l,
4470e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4480e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
449b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %gi",bs*a->j[k]+l,
4500e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4510e6d2581SBarry Smith             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
452b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4530ef38995SBarry Smith             }
45444cd7ae7SLois Curfman McInnes #else
4550ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
456b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
4570ef38995SBarry Smith             }
45844cd7ae7SLois Curfman McInnes #endif
45944cd7ae7SLois Curfman McInnes           }
46044cd7ae7SLois Curfman McInnes         }
461b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
46244cd7ae7SLois Curfman McInnes       }
46344cd7ae7SLois Curfman McInnes     }
464b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
4650ef38995SBarry Smith   } else {
466b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
467b6490206SBarry Smith     for (i=0; i<a->mbs; i++) {
468b6490206SBarry Smith       for (j=0; j<bs; j++) {
469b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
470b6490206SBarry Smith         for (k=a->i[i]; k<a->i[i+1]; k++) {
471b6490206SBarry Smith           for (l=0; l<bs; l++) {
472aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4730e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
474b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
4750e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4760e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
477b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
4780e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4790ef38995SBarry Smith             } else {
480b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
48188685aaeSLois Curfman McInnes             }
48288685aaeSLois Curfman McInnes #else
483b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
48488685aaeSLois Curfman McInnes #endif
4852593348eSBarry Smith           }
4862593348eSBarry Smith         }
487b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
4882593348eSBarry Smith       }
4892593348eSBarry Smith     }
490b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
491b6490206SBarry Smith   }
492b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
4933a40ed3dSBarry Smith   PetscFunctionReturn(0);
4942593348eSBarry Smith }
4952593348eSBarry Smith 
4964a2ae208SSatish Balay #undef __FUNCT__
4974a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom"
498b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
4993270192aSSatish Balay {
50077ed5343SBarry Smith   Mat          A = (Mat) Aa;
5013270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ*)A->data;
502b65db4caSBarry Smith   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2;
5030e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
5043f1db9ecSBarry Smith   MatScalar    *aa;
505b0a32e0cSBarry Smith   PetscViewer  viewer;
5063270192aSSatish Balay 
5073a40ed3dSBarry Smith   PetscFunctionBegin;
5083270192aSSatish Balay 
509b65db4caSBarry Smith   /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/
51077ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
51177ed5343SBarry Smith 
512b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
51377ed5343SBarry Smith 
5143270192aSSatish Balay   /* loop over matrix elements drawing boxes */
515b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
5163270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5173270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
518273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5193270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5203270192aSSatish Balay       aa = a->a + j*bs2;
5213270192aSSatish Balay       for (k=0; k<bs; k++) {
5223270192aSSatish Balay         for (l=0; l<bs; l++) {
5230e6d2581SBarry Smith           if (PetscRealPart(*aa++) >=  0.) continue;
524b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5253270192aSSatish Balay         }
5263270192aSSatish Balay       }
5273270192aSSatish Balay     }
5283270192aSSatish Balay   }
529b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
5303270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5313270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
532273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5333270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5343270192aSSatish Balay       aa = a->a + j*bs2;
5353270192aSSatish Balay       for (k=0; k<bs; k++) {
5363270192aSSatish Balay         for (l=0; l<bs; l++) {
5370e6d2581SBarry Smith           if (PetscRealPart(*aa++) != 0.) continue;
538b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5393270192aSSatish Balay         }
5403270192aSSatish Balay       }
5413270192aSSatish Balay     }
5423270192aSSatish Balay   }
5433270192aSSatish Balay 
544b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
5453270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5463270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
547273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5483270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5493270192aSSatish Balay       aa = a->a + j*bs2;
5503270192aSSatish Balay       for (k=0; k<bs; k++) {
5513270192aSSatish Balay         for (l=0; l<bs; l++) {
5520e6d2581SBarry Smith           if (PetscRealPart(*aa++) <= 0.) continue;
553b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5543270192aSSatish Balay         }
5553270192aSSatish Balay       }
5563270192aSSatish Balay     }
5573270192aSSatish Balay   }
55877ed5343SBarry Smith   PetscFunctionReturn(0);
55977ed5343SBarry Smith }
5603270192aSSatish Balay 
5614a2ae208SSatish Balay #undef __FUNCT__
5624a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw"
563b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer)
56477ed5343SBarry Smith {
5657dce120fSSatish Balay   int          ierr;
5660e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,w,h;
567b0a32e0cSBarry Smith   PetscDraw    draw;
56877ed5343SBarry Smith   PetscTruth   isnull;
5693270192aSSatish Balay 
57077ed5343SBarry Smith   PetscFunctionBegin;
57177ed5343SBarry Smith 
572b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
573b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
57477ed5343SBarry Smith 
57577ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
576273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
57777ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
578b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
579b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
58077ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5813a40ed3dSBarry Smith   PetscFunctionReturn(0);
5823270192aSSatish Balay }
5833270192aSSatish Balay 
5844a2ae208SSatish Balay #undef __FUNCT__
5854a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ"
586b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer)
5872593348eSBarry Smith {
58819bcc07fSBarry Smith   int        ierr;
5896831982aSBarry Smith   PetscTruth issocket,isascii,isbinary,isdraw;
5902593348eSBarry Smith 
5913a40ed3dSBarry Smith   PetscFunctionBegin;
592b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
593b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
594fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
595fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
5960f5bd95cSBarry Smith   if (issocket) {
59729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
5980f5bd95cSBarry Smith   } else if (isascii){
5993a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
6000f5bd95cSBarry Smith   } else if (isbinary) {
6013a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
6020f5bd95cSBarry Smith   } else if (isdraw) {
6033a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
6045cd90555SBarry Smith   } else {
60529bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
6062593348eSBarry Smith   }
6073a40ed3dSBarry Smith   PetscFunctionReturn(0);
6082593348eSBarry Smith }
609b6490206SBarry Smith 
610cd0e1443SSatish Balay 
6114a2ae208SSatish Balay #undef __FUNCT__
6124a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqBAIJ"
61387828ca2SBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v)
614cd0e1443SSatish Balay {
615cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6162d61bbb3SSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
6172d61bbb3SSatish Balay   int        *ai = a->i,*ailen = a->ilen;
6182d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6193f1db9ecSBarry Smith   MatScalar  *ap,*aa = a->a,zero = 0.0;
620cd0e1443SSatish Balay 
6213a40ed3dSBarry Smith   PetscFunctionBegin;
6222d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
623cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
62429bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
625273d9f13SBarry Smith     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
6262d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6272c3acbe9SBarry Smith     nrow = ailen[brow];
6282d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
62929bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
630273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
6312d61bbb3SSatish Balay       col  = in[l] ;
6322d61bbb3SSatish Balay       bcol = col/bs;
6332d61bbb3SSatish Balay       cidx = col%bs;
6342d61bbb3SSatish Balay       ridx = row%bs;
6352d61bbb3SSatish Balay       high = nrow;
6362d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6372d61bbb3SSatish Balay       while (high-low > 5) {
638cd0e1443SSatish Balay         t = (low+high)/2;
639cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
640cd0e1443SSatish Balay         else             low  = t;
641cd0e1443SSatish Balay       }
642cd0e1443SSatish Balay       for (i=low; i<high; i++) {
643cd0e1443SSatish Balay         if (rp[i] > bcol) break;
644cd0e1443SSatish Balay         if (rp[i] == bcol) {
6452d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6462d61bbb3SSatish Balay           goto finished;
647cd0e1443SSatish Balay         }
648cd0e1443SSatish Balay       }
6492d61bbb3SSatish Balay       *v++ = zero;
6502d61bbb3SSatish Balay       finished:;
651cd0e1443SSatish Balay     }
652cd0e1443SSatish Balay   }
6533a40ed3dSBarry Smith   PetscFunctionReturn(0);
654cd0e1443SSatish Balay }
655cd0e1443SSatish Balay 
65695d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
6574a2ae208SSatish Balay #undef __FUNCT__
6584a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
65987828ca2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,PetscScalar *v,InsertMode addv)
66095d5f7c2SBarry Smith {
661563b5814SBarry Smith   Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data;
662563b5814SBarry Smith   int         ierr,i,N = m*n*b->bs2;
663563b5814SBarry Smith   MatScalar   *vsingle;
66495d5f7c2SBarry Smith 
66595d5f7c2SBarry Smith   PetscFunctionBegin;
666563b5814SBarry Smith   if (N > b->setvalueslen) {
667563b5814SBarry Smith     if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);}
668b0a32e0cSBarry Smith     ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr);
669563b5814SBarry Smith     b->setvalueslen  = N;
670563b5814SBarry Smith   }
671563b5814SBarry Smith   vsingle = b->setvaluescopy;
67295d5f7c2SBarry Smith   for (i=0; i<N; i++) {
67395d5f7c2SBarry Smith     vsingle[i] = v[i];
67495d5f7c2SBarry Smith   }
67595d5f7c2SBarry Smith   ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr);
67695d5f7c2SBarry Smith   PetscFunctionReturn(0);
67795d5f7c2SBarry Smith }
67895d5f7c2SBarry Smith #endif
67995d5f7c2SBarry Smith 
6802d61bbb3SSatish Balay 
6814a2ae208SSatish Balay #undef __FUNCT__
6824a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
68395d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is)
68492c4ed94SBarry Smith {
68592c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6868a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
687273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
688549d3d68SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
689273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
69095d5f7c2SBarry Smith   MatScalar   *value = v,*ap,*aa = a->a,*bap;
69192c4ed94SBarry Smith 
6923a40ed3dSBarry Smith   PetscFunctionBegin;
6930e324ae4SSatish Balay   if (roworiented) {
6940e324ae4SSatish Balay     stepval = (n-1)*bs;
6950e324ae4SSatish Balay   } else {
6960e324ae4SSatish Balay     stepval = (m-1)*bs;
6970e324ae4SSatish Balay   }
69892c4ed94SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
69992c4ed94SBarry Smith     row  = im[k];
7005ef9f2a5SBarry Smith     if (row < 0) continue;
701aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
70229bbc08cSBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
70392c4ed94SBarry Smith #endif
70492c4ed94SBarry Smith     rp   = aj + ai[row];
70592c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
70692c4ed94SBarry Smith     rmax = imax[row];
70792c4ed94SBarry Smith     nrow = ailen[row];
70892c4ed94SBarry Smith     low  = 0;
70992c4ed94SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
7105ef9f2a5SBarry Smith       if (in[l] < 0) continue;
711aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
71229bbc08cSBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
71392c4ed94SBarry Smith #endif
71492c4ed94SBarry Smith       col = in[l];
71592c4ed94SBarry Smith       if (roworiented) {
7160e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
7170e324ae4SSatish Balay       } else {
7180e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
71992c4ed94SBarry Smith       }
72092c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
72192c4ed94SBarry Smith       while (high-low > 7) {
72292c4ed94SBarry Smith         t = (low+high)/2;
72392c4ed94SBarry Smith         if (rp[t] > col) high = t;
72492c4ed94SBarry Smith         else             low  = t;
72592c4ed94SBarry Smith       }
72692c4ed94SBarry Smith       for (i=low; i<high; i++) {
72792c4ed94SBarry Smith         if (rp[i] > col) break;
72892c4ed94SBarry Smith         if (rp[i] == col) {
7298a84c255SSatish Balay           bap  = ap +  bs2*i;
7300e324ae4SSatish Balay           if (roworiented) {
7318a84c255SSatish Balay             if (is == ADD_VALUES) {
732dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
733dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
7348a84c255SSatish Balay                   bap[jj] += *value++;
735dd9472c6SBarry Smith                 }
736dd9472c6SBarry Smith               }
7370e324ae4SSatish Balay             } else {
738dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
739dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
7400e324ae4SSatish Balay                   bap[jj] = *value++;
7418a84c255SSatish Balay                 }
742dd9472c6SBarry Smith               }
743dd9472c6SBarry Smith             }
7440e324ae4SSatish Balay           } else {
7450e324ae4SSatish Balay             if (is == ADD_VALUES) {
746dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
747dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
7480e324ae4SSatish Balay                   *bap++ += *value++;
749dd9472c6SBarry Smith                 }
750dd9472c6SBarry Smith               }
7510e324ae4SSatish Balay             } else {
752dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
753dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
7540e324ae4SSatish Balay                   *bap++  = *value++;
7550e324ae4SSatish Balay                 }
7568a84c255SSatish Balay               }
757dd9472c6SBarry Smith             }
758dd9472c6SBarry Smith           }
759f1241b54SBarry Smith           goto noinsert2;
76092c4ed94SBarry Smith         }
76192c4ed94SBarry Smith       }
76289280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
76329bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
76492c4ed94SBarry Smith       if (nrow >= rmax) {
76592c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
76692c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
7673f1db9ecSBarry Smith         MatScalar *new_a;
76892c4ed94SBarry Smith 
76929bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
77089280ab3SLois Curfman McInnes 
77192c4ed94SBarry Smith         /* malloc new storage space */
7723f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
773b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
77492c4ed94SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
77592c4ed94SBarry Smith         new_i   = new_j + new_nz;
77692c4ed94SBarry Smith 
77792c4ed94SBarry Smith         /* copy over old data into new slots */
77892c4ed94SBarry Smith         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
77992c4ed94SBarry Smith         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
780549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
78192c4ed94SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
782549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
783549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
784549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
785549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
78692c4ed94SBarry Smith         /* free up old matrix storage */
787606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
788606d414cSSatish Balay         if (!a->singlemalloc) {
789606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
790606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
791606d414cSSatish Balay         }
79292c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
793c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
79492c4ed94SBarry Smith 
79592c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
79692c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
797b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
79892c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
79992c4ed94SBarry Smith         a->reallocs++;
80092c4ed94SBarry Smith         a->nz++;
80192c4ed94SBarry Smith       }
80292c4ed94SBarry Smith       N = nrow++ - 1;
80392c4ed94SBarry Smith       /* shift up all the later entries in this row */
80492c4ed94SBarry Smith       for (ii=N; ii>=i; ii--) {
80592c4ed94SBarry Smith         rp[ii+1] = rp[ii];
806549d3d68SSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
80792c4ed94SBarry Smith       }
808549d3d68SSatish Balay       if (N >= i) {
809549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
810549d3d68SSatish Balay       }
81192c4ed94SBarry Smith       rp[i] = col;
8128a84c255SSatish Balay       bap   = ap +  bs2*i;
8130e324ae4SSatish Balay       if (roworiented) {
814dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
815dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs) {
8160e324ae4SSatish Balay             bap[jj] = *value++;
817dd9472c6SBarry Smith           }
818dd9472c6SBarry Smith         }
8190e324ae4SSatish Balay       } else {
820dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
821dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++) {
8220e324ae4SSatish Balay             *bap++  = *value++;
8230e324ae4SSatish Balay           }
824dd9472c6SBarry Smith         }
825dd9472c6SBarry Smith       }
826f1241b54SBarry Smith       noinsert2:;
82792c4ed94SBarry Smith       low = i;
82892c4ed94SBarry Smith     }
82992c4ed94SBarry Smith     ailen[row] = nrow;
83092c4ed94SBarry Smith   }
8313a40ed3dSBarry Smith   PetscFunctionReturn(0);
83292c4ed94SBarry Smith }
83392c4ed94SBarry Smith 
834*a56a16c8SHong Zhang EXTERN int MatUseDSCPACK_MPIBAIJ(Mat);
8354a2ae208SSatish Balay #undef __FUNCT__
8364a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ"
8378f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
838584200bdSSatish Balay {
839584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
840584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
841273d9f13SBarry Smith   int        m = A->m,*ip,N,*ailen = a->ilen;
842549d3d68SSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
8433f1db9ecSBarry Smith   MatScalar  *aa = a->a,*ap;
844*a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
845*a56a16c8SHong Zhang   PetscTruth   flag;
846*a56a16c8SHong Zhang #endif
847584200bdSSatish Balay 
8483a40ed3dSBarry Smith   PetscFunctionBegin;
8493a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
850584200bdSSatish Balay 
85143ee02c3SBarry Smith   if (m) rmax = ailen[0];
852584200bdSSatish Balay   for (i=1; i<mbs; i++) {
853584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
854584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
855d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
856584200bdSSatish Balay     if (fshift) {
857a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
858584200bdSSatish Balay       N = ailen[i];
859584200bdSSatish Balay       for (j=0; j<N; j++) {
860584200bdSSatish Balay         ip[j-fshift] = ip[j];
861549d3d68SSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
862584200bdSSatish Balay       }
863584200bdSSatish Balay     }
864584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
865584200bdSSatish Balay   }
866584200bdSSatish Balay   if (mbs) {
867584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
868584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
869584200bdSSatish Balay   }
870584200bdSSatish Balay   /* reset ilen and imax for each row */
871584200bdSSatish Balay   for (i=0; i<mbs; i++) {
872584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
873584200bdSSatish Balay   }
874a7c10996SSatish Balay   a->nz = ai[mbs];
875584200bdSSatish Balay 
876584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
877584200bdSSatish Balay   if (fshift && a->diag) {
878606d414cSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
879b424e231SHong Zhang     PetscLogObjectMemory(A,-(mbs+1)*sizeof(int));
880584200bdSSatish Balay     a->diag = 0;
881584200bdSSatish Balay   }
882bba1ac68SSatish Balay   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",m,A->n,a->bs,fshift*bs2,a->nz*bs2);
883bba1ac68SSatish Balay   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",a->reallocs);
884b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
885e2f3b5e9SSatish Balay   a->reallocs          = 0;
8860e6d2581SBarry Smith   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
887*a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
888*a56a16c8SHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
889*a56a16c8SHong Zhang   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(A);CHKERRQ(ierr); }
890*a56a16c8SHong Zhang #endif
891*a56a16c8SHong Zhang 
8924e220ebcSLois Curfman McInnes 
8933a40ed3dSBarry Smith   PetscFunctionReturn(0);
894584200bdSSatish Balay }
895584200bdSSatish Balay 
8965a838052SSatish Balay 
897bea157c4SSatish Balay 
898bea157c4SSatish Balay /*
899bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
900bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
901bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
902bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
903bea157c4SSatish Balay */
9044a2ae208SSatish Balay #undef __FUNCT__
9054a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks"
906bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
907d9b7c43dSSatish Balay {
908bea157c4SSatish Balay   int        i,j,k,row;
909bea157c4SSatish Balay   PetscTruth flg;
9103a40ed3dSBarry Smith 
911433994e6SBarry Smith   PetscFunctionBegin;
912bea157c4SSatish Balay   for (i=0,j=0; i<n; j++) {
913bea157c4SSatish Balay     row = idx[i];
914bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
915bea157c4SSatish Balay       sizes[j] = 1;
916bea157c4SSatish Balay       i++;
917e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
918bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
919bea157c4SSatish Balay       i++;
920bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
921bea157c4SSatish Balay       flg = PETSC_TRUE;
922bea157c4SSatish Balay       for (k=1; k<bs; k++) {
923bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
924bea157c4SSatish Balay           flg = PETSC_FALSE;
925bea157c4SSatish Balay           break;
926d9b7c43dSSatish Balay         }
927bea157c4SSatish Balay       }
928bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
929bea157c4SSatish Balay         sizes[j] = bs;
930bea157c4SSatish Balay         i+= bs;
931bea157c4SSatish Balay       } else {
932bea157c4SSatish Balay         sizes[j] = 1;
933bea157c4SSatish Balay         i++;
934bea157c4SSatish Balay       }
935bea157c4SSatish Balay     }
936bea157c4SSatish Balay   }
937bea157c4SSatish Balay   *bs_max = j;
9383a40ed3dSBarry Smith   PetscFunctionReturn(0);
939d9b7c43dSSatish Balay }
940d9b7c43dSSatish Balay 
9414a2ae208SSatish Balay #undef __FUNCT__
9424a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ"
94387828ca2SBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,PetscScalar *diag)
944d9b7c43dSSatish Balay {
945d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
946b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
947bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
94887828ca2SBarry Smith   PetscScalar zero = 0.0;
9493f1db9ecSBarry Smith   MatScalar   *aa;
950d9b7c43dSSatish Balay 
9513a40ed3dSBarry Smith   PetscFunctionBegin;
952d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
953b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr);
954d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
955d9b7c43dSSatish Balay 
956bea157c4SSatish Balay   /* allocate memory for rows,sizes */
957b0a32e0cSBarry Smith   ierr  = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
958bea157c4SSatish Balay   sizes = rows + is_n;
959bea157c4SSatish Balay 
960563b5814SBarry Smith   /* copy IS values to rows, and sort them */
961bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
962bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
963dffd3267SBarry Smith   if (baij->keepzeroedrows) {
964dffd3267SBarry Smith     for (i=0; i<is_n; i++) { sizes[i] = 1; }
965dffd3267SBarry Smith     bs_max = is_n;
966dffd3267SBarry Smith   } else {
967bea157c4SSatish Balay     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
968dffd3267SBarry Smith   }
969b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
970bea157c4SSatish Balay 
971bea157c4SSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
972bea157c4SSatish Balay     row   = rows[j];
973273d9f13SBarry Smith     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
974bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
975bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
976dffd3267SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
977bea157c4SSatish Balay       if (diag) {
978bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
979bea157c4SSatish Balay           baij->ilen[row/bs]       = 1;
980bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
981bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
982a07cd24cSSatish Balay         }
983563b5814SBarry Smith         /* Now insert all the diagonal values for this bs */
984bea157c4SSatish Balay         for (k=0; k<bs; k++) {
985bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
986bea157c4SSatish Balay         }
987bea157c4SSatish Balay       } else { /* (!diag) */
988bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
989bea157c4SSatish Balay       } /* end (!diag) */
990bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
991aa482453SBarry Smith #if defined (PETSC_USE_DEBUG)
99229bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
993bea157c4SSatish Balay #endif
994bea157c4SSatish Balay       for (k=0; k<count; k++) {
995d9b7c43dSSatish Balay         aa[0] =  zero;
996d9b7c43dSSatish Balay         aa    += bs;
997d9b7c43dSSatish Balay       }
998d9b7c43dSSatish Balay       if (diag) {
999f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
1000d9b7c43dSSatish Balay       }
1001d9b7c43dSSatish Balay     }
1002bea157c4SSatish Balay   }
1003bea157c4SSatish Balay 
1004606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
10059a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10063a40ed3dSBarry Smith   PetscFunctionReturn(0);
1007d9b7c43dSSatish Balay }
10081c351548SSatish Balay 
10094a2ae208SSatish Balay #undef __FUNCT__
10104a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ"
101187828ca2SBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
10122d61bbb3SSatish Balay {
10132d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
10142d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
1015273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
10162d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1017549d3d68SSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
1018273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
10193f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
10202d61bbb3SSatish Balay 
10212d61bbb3SSatish Balay   PetscFunctionBegin;
10222d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
10232d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
10245ef9f2a5SBarry Smith     if (row < 0) continue;
1025aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1026273d9f13SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
10272d61bbb3SSatish Balay #endif
10282d61bbb3SSatish Balay     rp   = aj + ai[brow];
10292d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
10302d61bbb3SSatish Balay     rmax = imax[brow];
10312d61bbb3SSatish Balay     nrow = ailen[brow];
10322d61bbb3SSatish Balay     low  = 0;
10332d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
10345ef9f2a5SBarry Smith       if (in[l] < 0) continue;
1035aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1036273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n);
10372d61bbb3SSatish Balay #endif
10382d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
10392d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
10402d61bbb3SSatish Balay       if (roworiented) {
10415ef9f2a5SBarry Smith         value = v[l + k*n];
10422d61bbb3SSatish Balay       } else {
10432d61bbb3SSatish Balay         value = v[k + l*m];
10442d61bbb3SSatish Balay       }
10452d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
10462d61bbb3SSatish Balay       while (high-low > 7) {
10472d61bbb3SSatish Balay         t = (low+high)/2;
10482d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
10492d61bbb3SSatish Balay         else              low  = t;
10502d61bbb3SSatish Balay       }
10512d61bbb3SSatish Balay       for (i=low; i<high; i++) {
10522d61bbb3SSatish Balay         if (rp[i] > bcol) break;
10532d61bbb3SSatish Balay         if (rp[i] == bcol) {
10542d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
10552d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
10562d61bbb3SSatish Balay           else                  *bap  = value;
10572d61bbb3SSatish Balay           goto noinsert1;
10582d61bbb3SSatish Balay         }
10592d61bbb3SSatish Balay       }
10602d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
106129bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
10622d61bbb3SSatish Balay       if (nrow >= rmax) {
10632d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
10642d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
10653f1db9ecSBarry Smith         MatScalar *new_a;
10662d61bbb3SSatish Balay 
106729bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
10682d61bbb3SSatish Balay 
10692d61bbb3SSatish Balay         /* Malloc new storage space */
10703f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
1071b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
10722d61bbb3SSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
10732d61bbb3SSatish Balay         new_i   = new_j + new_nz;
10742d61bbb3SSatish Balay 
10752d61bbb3SSatish Balay         /* copy over old data into new slots */
10762d61bbb3SSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
10772d61bbb3SSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1078549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
10792d61bbb3SSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1080549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1081549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1082549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1083549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
10842d61bbb3SSatish Balay         /* free up old matrix storage */
1085606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
1086606d414cSSatish Balay         if (!a->singlemalloc) {
1087606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
1088606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
1089606d414cSSatish Balay         }
10902d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1091c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
10922d61bbb3SSatish Balay 
10932d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10942d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1095b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
10962d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10972d61bbb3SSatish Balay         a->reallocs++;
10982d61bbb3SSatish Balay         a->nz++;
10992d61bbb3SSatish Balay       }
11002d61bbb3SSatish Balay       N = nrow++ - 1;
11012d61bbb3SSatish Balay       /* shift up all the later entries in this row */
11022d61bbb3SSatish Balay       for (ii=N; ii>=i; ii--) {
11032d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
1104549d3d68SSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
11052d61bbb3SSatish Balay       }
1106549d3d68SSatish Balay       if (N>=i) {
1107549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1108549d3d68SSatish Balay       }
11092d61bbb3SSatish Balay       rp[i]                      = bcol;
11102d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
11112d61bbb3SSatish Balay       noinsert1:;
11122d61bbb3SSatish Balay       low = i;
11132d61bbb3SSatish Balay     }
11142d61bbb3SSatish Balay     ailen[brow] = nrow;
11152d61bbb3SSatish Balay   }
11162d61bbb3SSatish Balay   PetscFunctionReturn(0);
11172d61bbb3SSatish Balay }
11182d61bbb3SSatish Balay 
11192d61bbb3SSatish Balay 
11204a2ae208SSatish Balay #undef __FUNCT__
11214a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ"
11225ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
11232d61bbb3SSatish Balay {
11242d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
11252d61bbb3SSatish Balay   Mat         outA;
11262d61bbb3SSatish Balay   int         ierr;
1127667159a5SBarry Smith   PetscTruth  row_identity,col_identity;
11282d61bbb3SSatish Balay 
11292d61bbb3SSatish Balay   PetscFunctionBegin;
113029bbc08cSBarry Smith   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
1131667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1132667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1133667159a5SBarry Smith   if (!row_identity || !col_identity) {
113429bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1135667159a5SBarry Smith   }
11362d61bbb3SSatish Balay 
11372d61bbb3SSatish Balay   outA          = inA;
11382d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
11392d61bbb3SSatish Balay 
11402d61bbb3SSatish Balay   if (!a->diag) {
1141c4992f7dSBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
11422d61bbb3SSatish Balay   }
1143cf242676SKris Buschelman 
1144c38d4ed2SBarry Smith   a->row        = row;
1145c38d4ed2SBarry Smith   a->col        = col;
1146c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1147c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1148c38d4ed2SBarry Smith 
1149c38d4ed2SBarry Smith   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
11504c49b128SBarry Smith   ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr);
1151b0a32e0cSBarry Smith   PetscLogObjectParent(inA,a->icol);
1152c38d4ed2SBarry Smith 
1153cf242676SKris Buschelman   /*
1154cf242676SKris Buschelman       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
1155cf242676SKris Buschelman       for ILU(0) factorization with natural ordering
1156cf242676SKris Buschelman   */
1157cf242676SKris Buschelman   if (a->bs < 8) {
1158cf242676SKris Buschelman     ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr);
1159cf242676SKris Buschelman   } else {
1160c38d4ed2SBarry Smith     if (!a->solve_work) {
116187828ca2SBarry Smith       ierr = PetscMalloc((inA->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
116287828ca2SBarry Smith       PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(PetscScalar));
1163c38d4ed2SBarry Smith     }
11642d61bbb3SSatish Balay   }
11652d61bbb3SSatish Balay 
1166667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1167667159a5SBarry Smith 
11682d61bbb3SSatish Balay   PetscFunctionReturn(0);
11692d61bbb3SSatish Balay }
11704a2ae208SSatish Balay #undef __FUNCT__
11714a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ"
1172ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1173ba4ca20aSSatish Balay {
1174c38d4ed2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1175ba4ca20aSSatish Balay   MPI_Comm          comm = A->comm;
1176d132466eSBarry Smith   int               ierr;
1177ba4ca20aSSatish Balay 
11783a40ed3dSBarry Smith   PetscFunctionBegin;
1179c38d4ed2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1180d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1181d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
11823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1183ba4ca20aSSatish Balay }
1184d9b7c43dSSatish Balay 
1185fb2e594dSBarry Smith EXTERN_C_BEGIN
11864a2ae208SSatish Balay #undef __FUNCT__
11874a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
118827a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
118927a8da17SBarry Smith {
119027a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
119114db4f74SSatish Balay   int         i,nz,nbs;
119227a8da17SBarry Smith 
119327a8da17SBarry Smith   PetscFunctionBegin;
119414db4f74SSatish Balay   nz  = baij->maxnz/baij->bs2;
119514db4f74SSatish Balay   nbs = baij->nbs;
119627a8da17SBarry Smith   for (i=0; i<nz; i++) {
119727a8da17SBarry Smith     baij->j[i] = indices[i];
119827a8da17SBarry Smith   }
119927a8da17SBarry Smith   baij->nz = nz;
120014db4f74SSatish Balay   for (i=0; i<nbs; i++) {
120127a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
120227a8da17SBarry Smith   }
120327a8da17SBarry Smith 
120427a8da17SBarry Smith   PetscFunctionReturn(0);
120527a8da17SBarry Smith }
1206fb2e594dSBarry Smith EXTERN_C_END
120727a8da17SBarry Smith 
12084a2ae208SSatish Balay #undef __FUNCT__
12094a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices"
121027a8da17SBarry Smith /*@
121127a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
121227a8da17SBarry Smith        in the matrix.
121327a8da17SBarry Smith 
121427a8da17SBarry Smith   Input Parameters:
121527a8da17SBarry Smith +  mat - the SeqBAIJ matrix
121627a8da17SBarry Smith -  indices - the column indices
121727a8da17SBarry Smith 
121815091d37SBarry Smith   Level: advanced
121915091d37SBarry Smith 
122027a8da17SBarry Smith   Notes:
122127a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
122227a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
122327a8da17SBarry Smith   of the MatSetValues() operation.
122427a8da17SBarry Smith 
122527a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
122627a8da17SBarry Smith   MatCreateSeqBAIJ().
122727a8da17SBarry Smith 
122827a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
122927a8da17SBarry Smith 
123027a8da17SBarry Smith @*/
123127a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
123227a8da17SBarry Smith {
123327a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
123427a8da17SBarry Smith 
123527a8da17SBarry Smith   PetscFunctionBegin;
123627a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1237c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
123827a8da17SBarry Smith   if (f) {
123927a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
124027a8da17SBarry Smith   } else {
124129bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
124227a8da17SBarry Smith   }
124327a8da17SBarry Smith   PetscFunctionReturn(0);
124427a8da17SBarry Smith }
124527a8da17SBarry Smith 
12464a2ae208SSatish Balay #undef __FUNCT__
12474a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ"
1248273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v)
1249273d9f13SBarry Smith {
1250273d9f13SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
1251273d9f13SBarry Smith   int          ierr,i,j,n,row,bs,*ai,*aj,mbs;
1252273d9f13SBarry Smith   PetscReal    atmp;
125387828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
1254273d9f13SBarry Smith   MatScalar    *aa;
1255273d9f13SBarry Smith   int          ncols,brow,krow,kcol;
1256273d9f13SBarry Smith 
1257273d9f13SBarry Smith   PetscFunctionBegin;
1258273d9f13SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1259273d9f13SBarry Smith   bs   = a->bs;
1260273d9f13SBarry Smith   aa   = a->a;
1261273d9f13SBarry Smith   ai   = a->i;
1262273d9f13SBarry Smith   aj   = a->j;
1263273d9f13SBarry Smith   mbs = a->mbs;
1264273d9f13SBarry Smith 
1265273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1266273d9f13SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1267273d9f13SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1268273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1269273d9f13SBarry Smith   for (i=0; i<mbs; i++) {
1270273d9f13SBarry Smith     ncols = ai[1] - ai[0]; ai++;
1271273d9f13SBarry Smith     brow  = bs*i;
1272273d9f13SBarry Smith     for (j=0; j<ncols; j++){
1273273d9f13SBarry Smith       /* bcol = bs*(*aj); */
1274273d9f13SBarry Smith       for (kcol=0; kcol<bs; kcol++){
1275273d9f13SBarry Smith         for (krow=0; krow<bs; krow++){
1276273d9f13SBarry Smith           atmp = PetscAbsScalar(*aa); aa++;
1277273d9f13SBarry Smith           row = brow + krow;    /* row index */
1278273d9f13SBarry Smith           /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */
1279273d9f13SBarry Smith           if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp;
1280273d9f13SBarry Smith         }
1281273d9f13SBarry Smith       }
1282273d9f13SBarry Smith       aj++;
1283273d9f13SBarry Smith     }
1284273d9f13SBarry Smith   }
1285273d9f13SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1286273d9f13SBarry Smith   PetscFunctionReturn(0);
1287273d9f13SBarry Smith }
1288273d9f13SBarry Smith 
12894a2ae208SSatish Balay #undef __FUNCT__
12904a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ"
1291273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A)
1292273d9f13SBarry Smith {
1293273d9f13SBarry Smith   int        ierr;
1294273d9f13SBarry Smith 
1295273d9f13SBarry Smith   PetscFunctionBegin;
1296273d9f13SBarry Smith   ierr =  MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1297273d9f13SBarry Smith   PetscFunctionReturn(0);
1298273d9f13SBarry Smith }
1299273d9f13SBarry Smith 
13004a2ae208SSatish Balay #undef __FUNCT__
13014a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ"
130287828ca2SBarry Smith int MatGetArray_SeqBAIJ(Mat A,PetscScalar **array)
1303f2a5309cSSatish Balay {
1304f2a5309cSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1305f2a5309cSSatish Balay   PetscFunctionBegin;
1306f2a5309cSSatish Balay   *array = a->a;
1307f2a5309cSSatish Balay   PetscFunctionReturn(0);
1308f2a5309cSSatish Balay }
1309f2a5309cSSatish Balay 
13104a2ae208SSatish Balay #undef __FUNCT__
13114a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ"
131287828ca2SBarry Smith int MatRestoreArray_SeqBAIJ(Mat A,PetscScalar **array)
1313f2a5309cSSatish Balay {
1314f2a5309cSSatish Balay   PetscFunctionBegin;
1315f2a5309cSSatish Balay   PetscFunctionReturn(0);
1316f2a5309cSSatish Balay }
1317f2a5309cSSatish Balay 
13182593348eSBarry Smith /* -------------------------------------------------------------------*/
1319cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1320cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1321cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1322cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1323cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
13247c922b88SBarry Smith        MatMultTranspose_SeqBAIJ,
13257c922b88SBarry Smith        MatMultTransposeAdd_SeqBAIJ,
1326cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1327cc2dc46cSBarry Smith        0,
1328cc2dc46cSBarry Smith        0,
1329cc2dc46cSBarry Smith        0,
1330cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1331cc2dc46cSBarry Smith        0,
1332b6490206SBarry Smith        0,
1333f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1334cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1335cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1336cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1337cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1338cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1339b6490206SBarry Smith        0,
1340cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1341cc2dc46cSBarry Smith        0,
1342cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1343cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1344cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1345cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1346cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1347cc2dc46cSBarry Smith        0,
1348cc2dc46cSBarry Smith        0,
1349273d9f13SBarry Smith        MatSetUpPreallocation_SeqBAIJ,
1350cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1351cc2dc46cSBarry Smith        0,
1352f2a5309cSSatish Balay        MatGetArray_SeqBAIJ,
1353f2a5309cSSatish Balay        MatRestoreArray_SeqBAIJ,
13542e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1355cc2dc46cSBarry Smith        0,
1356cc2dc46cSBarry Smith        0,
1357cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1358cc2dc46cSBarry Smith        0,
1359cc2dc46cSBarry Smith        0,
1360cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1361cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1362cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1363cc2dc46cSBarry Smith        0,
1364cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1365cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1366cc2dc46cSBarry Smith        0,
1367cc2dc46cSBarry Smith        0,
1368cc2dc46cSBarry Smith        0,
1369cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
13703b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
137192c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1372cc2dc46cSBarry Smith        0,
1373cc2dc46cSBarry Smith        0,
1374cc2dc46cSBarry Smith        0,
1375cc2dc46cSBarry Smith        0,
1376cc2dc46cSBarry Smith        0,
1377cc2dc46cSBarry Smith        0,
1378d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1379cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1380b9b97703SBarry Smith        MatDestroy_SeqBAIJ,
1381b9b97703SBarry Smith        MatView_SeqBAIJ,
13823a64cc32SBarry Smith        MatGetPetscMaps_Petsc,
1383273d9f13SBarry Smith        0,
1384273d9f13SBarry Smith        0,
1385273d9f13SBarry Smith        0,
1386273d9f13SBarry Smith        0,
1387273d9f13SBarry Smith        0,
1388273d9f13SBarry Smith        0,
1389273d9f13SBarry Smith        MatGetRowMax_SeqBAIJ,
1390273d9f13SBarry Smith        MatConvert_Basic};
13912593348eSBarry Smith 
13923e90b805SBarry Smith EXTERN_C_BEGIN
13934a2ae208SSatish Balay #undef __FUNCT__
13944a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ"
13953e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
13963e90b805SBarry Smith {
13973e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1398273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1399d132466eSBarry Smith   int         ierr;
14003e90b805SBarry Smith 
14013e90b805SBarry Smith   PetscFunctionBegin;
14023e90b805SBarry Smith   if (aij->nonew != 1) {
140329bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
14043e90b805SBarry Smith   }
14053e90b805SBarry Smith 
14063e90b805SBarry Smith   /* allocate space for values if not already there */
14073e90b805SBarry Smith   if (!aij->saved_values) {
140887828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
14093e90b805SBarry Smith   }
14103e90b805SBarry Smith 
14113e90b805SBarry Smith   /* copy values over */
141287828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
14133e90b805SBarry Smith   PetscFunctionReturn(0);
14143e90b805SBarry Smith }
14153e90b805SBarry Smith EXTERN_C_END
14163e90b805SBarry Smith 
14173e90b805SBarry Smith EXTERN_C_BEGIN
14184a2ae208SSatish Balay #undef __FUNCT__
14194a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ"
142032e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
14213e90b805SBarry Smith {
14223e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1423273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
14243e90b805SBarry Smith 
14253e90b805SBarry Smith   PetscFunctionBegin;
14263e90b805SBarry Smith   if (aij->nonew != 1) {
142729bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
14283e90b805SBarry Smith   }
14293e90b805SBarry Smith   if (!aij->saved_values) {
143029bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
14313e90b805SBarry Smith   }
14323e90b805SBarry Smith 
14333e90b805SBarry Smith   /* copy values over */
143487828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
14353e90b805SBarry Smith   PetscFunctionReturn(0);
14363e90b805SBarry Smith }
14373e90b805SBarry Smith EXTERN_C_END
14383e90b805SBarry Smith 
1439273d9f13SBarry Smith EXTERN_C_BEGIN
1440273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*);
1441273d9f13SBarry Smith EXTERN_C_END
1442273d9f13SBarry Smith 
1443273d9f13SBarry Smith EXTERN_C_BEGIN
14444a2ae208SSatish Balay #undef __FUNCT__
14454a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ"
1446273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B)
14472593348eSBarry Smith {
1448273d9f13SBarry Smith   int         ierr,size;
1449b6490206SBarry Smith   Mat_SeqBAIJ *b;
14503b2fbd54SBarry Smith 
14513a40ed3dSBarry Smith   PetscFunctionBegin;
1452273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
145329bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1454b6490206SBarry Smith 
1455273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1456273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1457b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr);
1458b0a32e0cSBarry Smith   B->data = (void*)b;
1459549d3d68SSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1460549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
14612593348eSBarry Smith   B->factor           = 0;
14622593348eSBarry Smith   B->lupivotthreshold = 1.0;
146390f02eecSBarry Smith   B->mapping          = 0;
14642593348eSBarry Smith   b->row              = 0;
14652593348eSBarry Smith   b->col              = 0;
1466e51c0b9cSSatish Balay   b->icol             = 0;
14672593348eSBarry Smith   b->reallocs         = 0;
14683e90b805SBarry Smith   b->saved_values     = 0;
1469cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
1470563b5814SBarry Smith   b->setvalueslen     = 0;
1471563b5814SBarry Smith   b->setvaluescopy    = PETSC_NULL;
1472563b5814SBarry Smith #endif
14738661488fSKris Buschelman   b->single_precision_solves = PETSC_FALSE;
14742593348eSBarry Smith 
14753a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
14763a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1477a5ae1ecdSBarry Smith 
1478c4992f7dSBarry Smith   b->sorted           = PETSC_FALSE;
1479c4992f7dSBarry Smith   b->roworiented      = PETSC_TRUE;
14802593348eSBarry Smith   b->nonew            = 0;
14812593348eSBarry Smith   b->diag             = 0;
14822593348eSBarry Smith   b->solve_work       = 0;
1483de6a44a3SBarry Smith   b->mult_work        = 0;
14842593348eSBarry Smith   b->spptr            = 0;
14850e6d2581SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1486883fce79SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
14874e220ebcSLois Curfman McInnes 
1488f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
14893e90b805SBarry Smith                                      "MatStoreValues_SeqBAIJ",
1490bc4b532fSSatish Balay                                       MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1491f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
14923e90b805SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
1493bc4b532fSSatish Balay                                       MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1494f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
149527a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
1496bc4b532fSSatish Balay                                       MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
1497273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_SeqBAIJ_SeqAIJ_C",
1498273d9f13SBarry Smith                                      "MatConvert_SeqBAIJ_SeqAIJ",
1499273d9f13SBarry Smith                                       MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr);
15003a40ed3dSBarry Smith   PetscFunctionReturn(0);
15012593348eSBarry Smith }
1502273d9f13SBarry Smith EXTERN_C_END
15032593348eSBarry Smith 
15044a2ae208SSatish Balay #undef __FUNCT__
15054a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ"
15062e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
15072593348eSBarry Smith {
15082593348eSBarry Smith   Mat         C;
1509b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data;
151027a8da17SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1511de6a44a3SBarry Smith 
15123a40ed3dSBarry Smith   PetscFunctionBegin;
151329bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
15142593348eSBarry Smith 
15152593348eSBarry Smith   *B = 0;
1516273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1517273d9f13SBarry Smith   ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr);
1518273d9f13SBarry Smith   c    = (Mat_SeqBAIJ*)C->data;
151944cd7ae7SLois Curfman McInnes 
1520b6490206SBarry Smith   c->bs         = a->bs;
15219df24120SSatish Balay   c->bs2        = a->bs2;
1522b6490206SBarry Smith   c->mbs        = a->mbs;
1523b6490206SBarry Smith   c->nbs        = a->nbs;
152435d8aa7fSBarry Smith   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
15252593348eSBarry Smith 
1526b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1527b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
1528b6490206SBarry Smith   for (i=0; i<mbs; i++) {
15292593348eSBarry Smith     c->imax[i] = a->imax[i];
15302593348eSBarry Smith     c->ilen[i] = a->ilen[i];
15312593348eSBarry Smith   }
15322593348eSBarry Smith 
15332593348eSBarry Smith   /* allocate the matrix space */
1534c4992f7dSBarry Smith   c->singlemalloc = PETSC_TRUE;
15353f1db9ecSBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1536b0a32e0cSBarry Smith   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
15377e67e3f9SSatish Balay   c->j = (int*)(c->a + nz*bs2);
1538de6a44a3SBarry Smith   c->i = c->j + nz;
1539549d3d68SSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1540b6490206SBarry Smith   if (mbs > 0) {
1541549d3d68SSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
15422e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1543549d3d68SSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
15442e8a6d31SBarry Smith     } else {
1545549d3d68SSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
15462593348eSBarry Smith     }
15472593348eSBarry Smith   }
15482593348eSBarry Smith 
1549b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
15502593348eSBarry Smith   c->sorted      = a->sorted;
15512593348eSBarry Smith   c->roworiented = a->roworiented;
15522593348eSBarry Smith   c->nonew       = a->nonew;
15532593348eSBarry Smith 
15542593348eSBarry Smith   if (a->diag) {
1555b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1556b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
1557b6490206SBarry Smith     for (i=0; i<mbs; i++) {
15582593348eSBarry Smith       c->diag[i] = a->diag[i];
15592593348eSBarry Smith     }
156098305bb5SBarry Smith   } else c->diag        = 0;
15612593348eSBarry Smith   c->nz                 = a->nz;
15622593348eSBarry Smith   c->maxnz              = a->maxnz;
15632593348eSBarry Smith   c->solve_work         = 0;
15642593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
15657fc0212eSBarry Smith   c->mult_work          = 0;
1566273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
1567273d9f13SBarry Smith   C->assembled          = PETSC_TRUE;
15682593348eSBarry Smith   *B = C;
1569b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
15703a40ed3dSBarry Smith   PetscFunctionReturn(0);
15712593348eSBarry Smith }
15722593348eSBarry Smith 
1573273d9f13SBarry Smith EXTERN_C_BEGIN
15744a2ae208SSatish Balay #undef __FUNCT__
15754a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ"
1576b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A)
15772593348eSBarry Smith {
1578b6490206SBarry Smith   Mat_SeqBAIJ  *a;
15792593348eSBarry Smith   Mat          B;
1580f1af5d2fSBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1581b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
158235aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1583a7c10996SSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
158487828ca2SBarry Smith   PetscScalar  *aa;
158519bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
15862593348eSBarry Smith 
15873a40ed3dSBarry Smith   PetscFunctionBegin;
1588b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1589de6a44a3SBarry Smith   bs2  = bs*bs;
1590b6490206SBarry Smith 
1591d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
159229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1593b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
15940752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1595552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
15962593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
15972593348eSBarry Smith 
1598d64ed03dSBarry Smith   if (header[3] < 0) {
159929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ");
1600d64ed03dSBarry Smith   }
1601d64ed03dSBarry Smith 
160229bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
160335aab85fSBarry Smith 
160435aab85fSBarry Smith   /*
160535aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
160635aab85fSBarry Smith     divisible by the blocksize
160735aab85fSBarry Smith   */
1608b6490206SBarry Smith   mbs        = M/bs;
160935aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
161035aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
161135aab85fSBarry Smith   else                  mbs++;
161235aab85fSBarry Smith   if (extra_rows) {
1613b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
161435aab85fSBarry Smith   }
1615b6490206SBarry Smith 
16162593348eSBarry Smith   /* read in row lengths */
1617b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
16180752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
161935aab85fSBarry Smith   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
16202593348eSBarry Smith 
1621b6490206SBarry Smith   /* read in column indices */
1622b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
16230752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
162435aab85fSBarry Smith   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
1625b6490206SBarry Smith 
1626b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1627b0a32e0cSBarry Smith   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
1628549d3d68SSatish Balay   ierr     = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
1629b0a32e0cSBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
1630549d3d68SSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
163135aab85fSBarry Smith   masked   = mask + mbs;
1632b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1633b6490206SBarry Smith   for (i=0; i<mbs; i++) {
163435aab85fSBarry Smith     nmask = 0;
1635b6490206SBarry Smith     for (j=0; j<bs; j++) {
1636b6490206SBarry Smith       kmax = rowlengths[rowcount];
1637b6490206SBarry Smith       for (k=0; k<kmax; k++) {
163835aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
163935aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1640b6490206SBarry Smith       }
1641b6490206SBarry Smith       rowcount++;
1642b6490206SBarry Smith     }
164335aab85fSBarry Smith     browlengths[i] += nmask;
164435aab85fSBarry Smith     /* zero out the mask elements we set */
164535aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1646b6490206SBarry Smith   }
1647b6490206SBarry Smith 
16482593348eSBarry Smith   /* create our matrix */
16493eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
16502593348eSBarry Smith   B = *A;
1651b6490206SBarry Smith   a = (Mat_SeqBAIJ*)B->data;
16522593348eSBarry Smith 
16532593348eSBarry Smith   /* set matrix "i" values */
1654de6a44a3SBarry Smith   a->i[0] = 0;
1655b6490206SBarry Smith   for (i=1; i<= mbs; i++) {
1656b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1657b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
16582593348eSBarry Smith   }
1659b6490206SBarry Smith   a->nz         = 0;
1660b6490206SBarry Smith   for (i=0; i<mbs; i++) a->nz += browlengths[i];
16612593348eSBarry Smith 
1662b6490206SBarry Smith   /* read in nonzero values */
166387828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
16640752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
166535aab85fSBarry Smith   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
1666b6490206SBarry Smith 
1667b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1668b6490206SBarry Smith   nzcount = 0; jcount = 0;
1669b6490206SBarry Smith   for (i=0; i<mbs; i++) {
1670b6490206SBarry Smith     nzcountb = nzcount;
167135aab85fSBarry Smith     nmask    = 0;
1672b6490206SBarry Smith     for (j=0; j<bs; j++) {
1673b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1674b6490206SBarry Smith       for (k=0; k<kmax; k++) {
167535aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
167635aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1677b6490206SBarry Smith       }
1678b6490206SBarry Smith     }
1679de6a44a3SBarry Smith     /* sort the masked values */
1680433994e6SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1681de6a44a3SBarry Smith 
1682b6490206SBarry Smith     /* set "j" values into matrix */
1683b6490206SBarry Smith     maskcount = 1;
168435aab85fSBarry Smith     for (j=0; j<nmask; j++) {
168535aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1686de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1687b6490206SBarry Smith     }
1688b6490206SBarry Smith     /* set "a" values into matrix */
1689de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1690b6490206SBarry Smith     for (j=0; j<bs; j++) {
1691b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1692b6490206SBarry Smith       for (k=0; k<kmax; k++) {
1693de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1694de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1695de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1696de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1697375fe846SBarry Smith         a->a[idx] = (MatScalar)aa[nzcountb++];
1698b6490206SBarry Smith       }
1699b6490206SBarry Smith     }
170035aab85fSBarry Smith     /* zero out the mask elements we set */
170135aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1702b6490206SBarry Smith   }
170329bbc08cSBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
1704b6490206SBarry Smith 
1705606d414cSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1706606d414cSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1707606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
1708606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
1709606d414cSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
1710b6490206SBarry Smith 
1711b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1712de6a44a3SBarry Smith 
17139c01be13SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
17143a40ed3dSBarry Smith   PetscFunctionReturn(0);
17152593348eSBarry Smith }
1716273d9f13SBarry Smith EXTERN_C_END
17172593348eSBarry Smith 
17184a2ae208SSatish Balay #undef __FUNCT__
17194a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ"
1720273d9f13SBarry Smith /*@C
1721273d9f13SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
1722273d9f13SBarry Smith    compressed row) format.  For good matrix assembly performance the
1723273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1724273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1725273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
17262593348eSBarry Smith 
1727273d9f13SBarry Smith    Collective on MPI_Comm
1728273d9f13SBarry Smith 
1729273d9f13SBarry Smith    Input Parameters:
1730273d9f13SBarry Smith +  comm - MPI communicator, set to PETSC_COMM_SELF
1731273d9f13SBarry Smith .  bs - size of block
1732273d9f13SBarry Smith .  m - number of rows
1733273d9f13SBarry Smith .  n - number of columns
173435d8aa7fSBarry Smith .  nz - number of nonzero blocks  per block row (same for all rows)
173535d8aa7fSBarry Smith -  nnz - array containing the number of nonzero blocks in the various block rows
1736273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1737273d9f13SBarry Smith 
1738273d9f13SBarry Smith    Output Parameter:
1739273d9f13SBarry Smith .  A - the matrix
1740273d9f13SBarry Smith 
1741273d9f13SBarry Smith    Options Database Keys:
1742273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1743273d9f13SBarry Smith                      block calculations (much slower)
1744273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1745273d9f13SBarry Smith 
1746273d9f13SBarry Smith    Level: intermediate
1747273d9f13SBarry Smith 
1748273d9f13SBarry Smith    Notes:
174935d8aa7fSBarry Smith    A nonzero block is any block that as 1 or more nonzeros in it
175035d8aa7fSBarry Smith 
1751273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1752273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1753273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1754273d9f13SBarry Smith 
1755273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1756273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1757273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1758273d9f13SBarry Smith    matrices.
1759273d9f13SBarry Smith 
1760273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1761273d9f13SBarry Smith @*/
1762273d9f13SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1763273d9f13SBarry Smith {
1764273d9f13SBarry Smith   int ierr;
1765273d9f13SBarry Smith 
1766273d9f13SBarry Smith   PetscFunctionBegin;
1767273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1768273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr);
1769273d9f13SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1770273d9f13SBarry Smith   PetscFunctionReturn(0);
1771273d9f13SBarry Smith }
1772273d9f13SBarry Smith 
17734a2ae208SSatish Balay #undef __FUNCT__
17744a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation"
1775273d9f13SBarry Smith /*@C
1776273d9f13SBarry Smith    MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros
1777273d9f13SBarry Smith    per row in the matrix. For good matrix assembly performance the
1778273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1779273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1780273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
1781273d9f13SBarry Smith 
1782273d9f13SBarry Smith    Collective on MPI_Comm
1783273d9f13SBarry Smith 
1784273d9f13SBarry Smith    Input Parameters:
1785273d9f13SBarry Smith +  A - the matrix
1786273d9f13SBarry Smith .  bs - size of block
1787273d9f13SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1788273d9f13SBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
1789273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1790273d9f13SBarry Smith 
1791273d9f13SBarry Smith    Options Database Keys:
1792273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1793273d9f13SBarry Smith                      block calculations (much slower)
1794273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1795273d9f13SBarry Smith 
1796273d9f13SBarry Smith    Level: intermediate
1797273d9f13SBarry Smith 
1798273d9f13SBarry Smith    Notes:
1799273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1800273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1801273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1802273d9f13SBarry Smith 
1803273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1804273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1805273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1806273d9f13SBarry Smith    matrices.
1807273d9f13SBarry Smith 
1808273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1809273d9f13SBarry Smith @*/
1810273d9f13SBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
1811273d9f13SBarry Smith {
1812273d9f13SBarry Smith   Mat_SeqBAIJ *b;
1813273d9f13SBarry Smith   int         i,len,ierr,mbs,nbs,bs2,newbs = bs;
1814273d9f13SBarry Smith   PetscTruth  flg;
1815273d9f13SBarry Smith 
1816273d9f13SBarry Smith   PetscFunctionBegin;
1817273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flg);CHKERRQ(ierr);
1818273d9f13SBarry Smith   if (!flg) PetscFunctionReturn(0);
1819273d9f13SBarry Smith 
1820273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1821b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr);
1822273d9f13SBarry Smith   if (nnz && newbs != bs) {
1823273d9f13SBarry Smith     SETERRQ(1,"Cannot change blocksize from command line if setting nnz");
1824273d9f13SBarry Smith   }
1825273d9f13SBarry Smith   bs   = newbs;
1826273d9f13SBarry Smith 
1827273d9f13SBarry Smith   mbs  = B->m/bs;
1828273d9f13SBarry Smith   nbs  = B->n/bs;
1829273d9f13SBarry Smith   bs2  = bs*bs;
1830273d9f13SBarry Smith 
1831273d9f13SBarry Smith   if (mbs*bs!=B->m || nbs*bs!=B->n) {
183235d8aa7fSBarry Smith     SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs);
1833273d9f13SBarry Smith   }
1834273d9f13SBarry Smith 
1835435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1836435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
1837273d9f13SBarry Smith   if (nnz) {
1838273d9f13SBarry Smith     for (i=0; i<mbs; i++) {
1839273d9f13SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
18403a7fca6bSBarry 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);
1841273d9f13SBarry Smith     }
1842273d9f13SBarry Smith   }
1843273d9f13SBarry Smith 
1844273d9f13SBarry Smith   b       = (Mat_SeqBAIJ*)B->data;
1845b0a32e0cSBarry Smith   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
184654138f6bSKris Buschelman   B->ops->solve               = MatSolve_SeqBAIJ_Update;
184754138f6bSKris Buschelman   B->ops->solvetranspose      = MatSolveTranspose_SeqBAIJ_Update;
1848273d9f13SBarry Smith   if (!flg) {
1849273d9f13SBarry Smith     switch (bs) {
1850273d9f13SBarry Smith     case 1:
1851273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1852273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1853273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
1854273d9f13SBarry Smith       break;
1855273d9f13SBarry Smith     case 2:
1856273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1857273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1858273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
1859273d9f13SBarry Smith       break;
1860273d9f13SBarry Smith     case 3:
1861273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1862273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1863273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
1864273d9f13SBarry Smith       break;
1865273d9f13SBarry Smith     case 4:
1866273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1867273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1868273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
1869273d9f13SBarry Smith       break;
1870273d9f13SBarry Smith     case 5:
1871273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1872273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1873273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
1874273d9f13SBarry Smith       break;
1875273d9f13SBarry Smith     case 6:
1876273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
1877273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
1878273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
1879273d9f13SBarry Smith       break;
1880273d9f13SBarry Smith     case 7:
188154138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
1882273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1883273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
1884273d9f13SBarry Smith       break;
1885273d9f13SBarry Smith     default:
188654138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
1887273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_N;
1888273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_N;
1889273d9f13SBarry Smith       break;
1890273d9f13SBarry Smith     }
1891273d9f13SBarry Smith   }
1892273d9f13SBarry Smith   b->bs      = bs;
1893273d9f13SBarry Smith   b->mbs     = mbs;
1894273d9f13SBarry Smith   b->nbs     = nbs;
1895b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
1896273d9f13SBarry Smith   if (!nnz) {
1897435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1898273d9f13SBarry Smith     else if (nz <= 0)        nz = 1;
1899273d9f13SBarry Smith     for (i=0; i<mbs; i++) b->imax[i] = nz;
1900273d9f13SBarry Smith     nz = nz*mbs;
1901273d9f13SBarry Smith   } else {
1902273d9f13SBarry Smith     nz = 0;
1903273d9f13SBarry Smith     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1904273d9f13SBarry Smith   }
1905273d9f13SBarry Smith 
1906273d9f13SBarry Smith   /* allocate the matrix space */
1907273d9f13SBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
1908b0a32e0cSBarry Smith   ierr  = PetscMalloc(len,&b->a);CHKERRQ(ierr);
1909273d9f13SBarry Smith   ierr  = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1910273d9f13SBarry Smith   b->j  = (int*)(b->a + nz*bs2);
1911273d9f13SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
1912273d9f13SBarry Smith   b->i  = b->j + nz;
1913273d9f13SBarry Smith   b->singlemalloc = PETSC_TRUE;
1914273d9f13SBarry Smith 
1915273d9f13SBarry Smith   b->i[0] = 0;
1916273d9f13SBarry Smith   for (i=1; i<mbs+1; i++) {
1917273d9f13SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
1918273d9f13SBarry Smith   }
1919273d9f13SBarry Smith 
1920273d9f13SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
192116d6aa21SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1922b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1923273d9f13SBarry Smith   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
1924273d9f13SBarry Smith 
1925273d9f13SBarry Smith   b->bs               = bs;
1926273d9f13SBarry Smith   b->bs2              = bs2;
1927273d9f13SBarry Smith   b->mbs              = mbs;
1928273d9f13SBarry Smith   b->nz               = 0;
1929273d9f13SBarry Smith   b->maxnz            = nz*bs2;
1930273d9f13SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1931273d9f13SBarry Smith   PetscFunctionReturn(0);
1932273d9f13SBarry Smith }
19332593348eSBarry Smith 
1934